Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGenerator.cxx
Go to the documentation of this file.
1// @(#)root/eg:$Id$
2// Author: Ola Nordmann 21/09/95
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGenerator
13 \ingroup eg
14
15The interface to various event generators
16
17Is an base class, that defines the interface of ROOT to various
18event generators. Every event generator should inherit from
19TGenerator or its subclasses.
20
21Derived class can overload the member function GenerateEvent
22to do the actual event generation (e.g., call PYEVNT or similar).
23
24The derived class should overload the member function
25ImportParticles (both types) to read the internal storage of the
26generated event into either the internal TObjArray or the passed
27TClonesArray of TParticles.
28
29If the generator code stores event data in the /HEPEVT/ common block
30Then the default implementation of ImportParticles should suffice.
31The common block /HEPEVT/ is structed like
32
33\verbatim
34 // C
35 typedef struct {
36 Int_t nevhep; // Event number
37 Int_t nhep; // # of particles
38 Int_t isthep[4000]; // Status flag of i'th particle
39 Int_t idhep[4000]; // PDG # of particle
40 Int_t jmohep[4000][2]; // 1st & 2nd mother particle #
41 Int_t jdahep[4000][2]; // 1st & 2nd daughter particle #
42 Double_t phep[4000][5]; // 4-momentum and 1 word
43 Double_t vhep[4000][4]; // 4-position of production
44 } HEPEVT_DEF;
45
46
47 C Fortran
48 COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(4000),IDHEP(4000),
49 + JMOHEP(2,4000),JDAHEP(2,4000),PHEP(5,4000),VHEP(4,4000)
50 INTEGER NEVHEP,NHEP,ISTHEP,IDHEP,JMOHEP,JDAHEP
51 DOUBLE PRECISION PHEP,VHEP
52\endverbatim
53
54The generic member functions SetParameter and GetParameter can be
55overloaded to set and get parameters of the event generator.
56
57Note, if the derived class interfaces a (set of) Fortran common
58blocks (like TPythia, TVenus does), one better make the derived
59class a singleton. That is, something like
60
61\verbatim
62 class MyGenerator : public TGenerator
63 {
64 public:
65 static MyGenerator* Instance()
66 {
67 if (!fgInstance) fgInstance = new MyGenerator;
68 return fgInstance;
69 }
70 void GenerateEvent() { ... }
71 void ImportParticles(TClonesArray* a, Option_t opt="") {...}
72 Int_t ImportParticles(Option_t opt="") { ... }
73 Int_t SetParameter(const char* name, Double_t val) { ... }
74 Double_t GetParameter(const char* name) { ... }
75 virtual ~MyGenerator() { ... }
76 protected:
77 MyGenerator() { ... }
78 MyGenerator(const MyGenerator& o) { ... }
79 MyGenerator& operator=(const MyGenerator& o) { ... }
80 static MyGenerator* fgInstance;
81 ClassDefOverride(MyGenerator,0);
82 };
83\endverbatim
84
85Having multiple objects accessing the same common blocks is not
86safe.
87
88Concrete TGenerator classes can be loaded in scripts and subseqent-
89ly used in compiled code:
90
91\verbatim
92 // MyRun.h
93 class MyRun : public TObject
94 {
95 public:
96 static MyRun* Instance() { ... }
97 void SetGenerator(TGenerator* g) { fGenerator = g; }
98 void Run(Int_t n, Option_t* option="")
99 {
100 TFile* file = TFile::Open("file.root","RECREATE");
101 TTree* tree = new TTree("T","T");
102 TClonesArray* p = new TClonesArray("TParticles");
103 tree->Branch("particles", &p);
104 for (Int_t event = 0; event < n; event++) {
105 fGenerator->GenerateEvent();
106 fGenerator->ImportParticles(p,option);
107 tree->Fill();
108 }
109 file->Write();
110 file->Close();
111 }
112 ...
113 protected:
114 TGenerator* fGenerator;
115 ClassDefOverride(MyRun,0);
116 };
117
118 // Config.C
119 void Config()
120 {
121 MyRun* run = MyRun::Instance();
122 run->SetGenerator(MyGenerator::Instance());
123 }
124
125 // main.cxx
126 int
127 main(int argc, char** argv)
128 {
129 TApplication app("", 0, 0);
130 gSystem->ProcessLine(".x Config.C");
131 MyRun::Instance()->Run(10);
132 return 0;
133 }
134\endverbatim
135
136This is especially useful for example with TVirtualMC or similar.
137*/
138
139#include "TROOT.h"
140#include "TGenerator.h"
141#include "TDatabasePDG.h"
142#include "TParticlePDG.h"
143#include "TParticle.h"
144#include "TObjArray.h"
145#include "Hepevt.h"
146#include "TVirtualPad.h"
147#include "TView.h"
148#include "TText.h"
149#include "TPaveText.h"
150#include "TClonesArray.h"
151#include "strlcpy.h"
152#include "snprintf.h"
153
154#include <iostream>
155
156
158
159////////////////////////////////////////////////////////////////////////////////
160/// Event generator default constructor
161///
162
163TGenerator::TGenerator(const char *name,const char *title): TNamed(name,title)
164{
165 // Initialize particles table
167 //TDatabasePDG *pdg = TDatabasePDG::Instance();
168 //if (!pdg->ParticleList()) pdg->Init();
169
170 fPtCut = 0;
172 fParticles = new TObjArray(10000);
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Event generator default destructor
177///
178
180{
181 //do nothing
182 if (fParticles) {
184 delete fParticles;
185 fParticles = 0;
186 }
187}
188
189////////////////////////////////////////////////////////////////////////////////
190/// must be implemented in concrete class (see eg TPythia6)
191
193{
194}
195
196////////////////////////////////////////////////////////////////////////////////
197///
198/// It reads the /HEPEVT/ common block which has been filled by the
199/// GenerateEvent method. If the event generator does not use the
200/// HEPEVT common block, This routine has to be overloaded by the
201/// subclasses.
202///
203/// The default action is to store only the stable particles (ISTHEP =
204/// 1) This can be demanded explicitly by setting the option = "Final"
205/// If the option = "All", all the particles are stored.
206///
207
209{
210 fParticles->Clear();
211 Int_t numpart = HEPEVT.nhep;
212 if (!strcmp(option,"") || !strcmp(option,"Final")) {
213 for (Int_t i = 0; i<numpart; i++) {
214 if (HEPEVT.isthep[i] == 1) {
215//
216// Use the common block values for the TParticle constructor
217//
218 TParticle *p = new TParticle(
219 HEPEVT.idhep[i],
220 HEPEVT.isthep[i],
221 HEPEVT.jmohep[i][0]-1,
222 HEPEVT.jmohep[i][1]-1,
223 HEPEVT.jdahep[i][0]-1,
224 HEPEVT.jdahep[i][1]-1,
225 HEPEVT.phep[i][0],
226 HEPEVT.phep[i][1],
227 HEPEVT.phep[i][2],
228 HEPEVT.phep[i][3],
229 HEPEVT.vhep[i][0],
230 HEPEVT.vhep[i][1],
231 HEPEVT.vhep[i][2],
232 HEPEVT.vhep[i][3]);
233 fParticles->Add(p);
234 }
235 }
236 } else if (!strcmp(option,"All")) {
237 for (Int_t i = 0; i<numpart; i++) {
238 TParticle *p = new TParticle(
239 HEPEVT.idhep[i],
240 HEPEVT.isthep[i],
241 HEPEVT.jmohep[i][0]-1,
242 HEPEVT.jmohep[i][1]-1,
243 HEPEVT.jdahep[i][0]-1,
244 HEPEVT.jdahep[i][1]-1,
245 HEPEVT.phep[i][0],
246 HEPEVT.phep[i][1],
247 HEPEVT.phep[i][2],
248 HEPEVT.phep[i][3],
249 HEPEVT.vhep[i][0],
250 HEPEVT.vhep[i][1],
251 HEPEVT.vhep[i][2],
252 HEPEVT.vhep[i][3]);
253 fParticles->Add(p);
254 }
255 }
256 return fParticles;
257}
258
259////////////////////////////////////////////////////////////////////////////////
260///
261/// It reads the /HEPEVT/ common block which has been filled by the
262/// GenerateEvent method. If the event generator does not use the
263/// HEPEVT common block, This routine has to be overloaded by the
264/// subclasses.
265///
266/// The function loops on the generated particles and store them in
267/// the TClonesArray pointed by the argument particles. The default
268/// action is to store only the stable particles (ISTHEP = 1) This can
269/// be demanded explicitly by setting the option = "Final" If the
270/// option = "All", all the particles are stored.
271///
272
274{
275 if (particles == 0) return 0;
276 TClonesArray &clonesParticles = *particles;
277 clonesParticles.Clear();
278 Int_t numpart = HEPEVT.nhep;
279 if (!strcmp(option,"") || !strcmp(option,"Final")) {
280 for (Int_t i = 0; i<numpart; i++) {
281 if (HEPEVT.isthep[i] == 1) {
282//
283// Use the common block values for the TParticle constructor
284//
285 new(clonesParticles[i]) TParticle(
286 HEPEVT.idhep[i],
287 HEPEVT.isthep[i],
288 HEPEVT.jmohep[i][0]-1,
289 HEPEVT.jmohep[i][1]-1,
290 HEPEVT.jdahep[i][0]-1,
291 HEPEVT.jdahep[i][1]-1,
292 HEPEVT.phep[i][0],
293 HEPEVT.phep[i][1],
294 HEPEVT.phep[i][2],
295 HEPEVT.phep[i][3],
296 HEPEVT.vhep[i][0],
297 HEPEVT.vhep[i][1],
298 HEPEVT.vhep[i][2],
299 HEPEVT.vhep[i][3]);
300 }
301 }
302 } else if (!strcmp(option,"All")) {
303 for (Int_t i = 0; i<numpart; i++) {
304 new(clonesParticles[i]) TParticle(
305 HEPEVT.idhep[i],
306 HEPEVT.isthep[i],
307 HEPEVT.jmohep[i][0]-1,
308 HEPEVT.jmohep[i][1]-1,
309 HEPEVT.jdahep[i][0]-1,
310 HEPEVT.jdahep[i][1]-1,
311 HEPEVT.phep[i][0],
312 HEPEVT.phep[i][1],
313 HEPEVT.phep[i][2],
314 HEPEVT.phep[i][3],
315 HEPEVT.vhep[i][0],
316 HEPEVT.vhep[i][1],
317 HEPEVT.vhep[i][2],
318 HEPEVT.vhep[i][3]);
319 }
320 }
321 return numpart;
322}
323
324////////////////////////////////////////////////////////////////////////////////
325///browse generator
326
328{
329 Draw();
330 gPad->Update();
331}
332
333////////////////////////////////////////////////////////////////////////////////
334/// Compute distance from point px,py to objects in event
335///
336
338{
339 const Int_t big = 9999;
340 const Int_t inview = 0;
341 Int_t dist = big;
342 if (px > 50 && py > 50) dist = inview;
343 return dist;
344}
345
346////////////////////////////////////////////////////////////////////////////////
347///
348/// Insert one event in the pad list
349///
350
352{
353 // Create a default canvas if a canvas does not exist
354 if (!gPad) {
355 gROOT->MakeDefCanvas();
356 if (gPad->GetVirtCanvas())
357 gPad->GetVirtCanvas()->SetFillColor(13);
358 }
359
360 static Float_t rbox = 1000;
361 Float_t rmin[3],rmax[3];
362 TView *view = gPad->GetView();
363 if (!strstr(option,"same")) {
364 if (view) { view->GetRange(rmin,rmax); rbox = rmax[2];}
365 gPad->Clear();
366 }
367
369
370 view = gPad->GetView();
371 // compute 3D view
372 if (view) {
373 view->GetRange(rmin,rmax);
374 rbox = rmax[2];
375 } else {
376 view = TView::CreateView(1,0,0);
377 if (view) view->SetRange(-rbox,-rbox,-rbox, rbox,rbox,rbox );
378 }
379 const Int_t kColorProton = 4;
380 const Int_t kColorNeutron = 5;
381 const Int_t kColorAntiProton= 3;
382 const Int_t kColorPionPlus = 6;
383 const Int_t kColorPionMinus = 2;
384 const Int_t kColorKaons = 7;
385 const Int_t kColorElectrons = 0;
386 const Int_t kColorGamma = 18;
387
388 Int_t nProtons = 0;
389 Int_t nNeutrons = 0;
390 Int_t nAntiProtons= 0;
391 Int_t nPionPlus = 0;
392 Int_t nPionMinus = 0;
393 Int_t nKaons = 0;
394 Int_t nElectrons = 0;
395 Int_t nGammas = 0;
396
397 Int_t ntracks = fParticles->GetEntriesFast();
398 Int_t i,lwidth,color,lstyle;
399 TParticlePDG *ap;
400 TParticle *p;
401 const char *name;
402 Double_t etot,vx,vy,vz;
403 Int_t ninvol = 0;
404 for (i=0;i<ntracks;i++) {
406 if(!p) continue;
407 ap = (TParticlePDG*)p->GetPDG();
408 vx = p->Vx();
409 vy = p->Vy();
410 vz = p->Vz();
411 if (vx*vx+vy*vy+vz*vz > rbox*rbox) continue;
412 Float_t pt = p->Pt();
413 if (pt < fPtCut) continue;
414 etot = p->Energy();
415 if (etot > 0.1) lwidth = Int_t(6*TMath::Log10(etot));
416 else lwidth = 1;
417 if (lwidth < 1) lwidth = 1;
418 lstyle = 1;
419 color = 0;
420 name = ap->GetName();
421 if (!strcmp(name,"n")) { if (!fShowNeutrons) continue;
422 color = kColorNeutron; nNeutrons++;}
423 if (!strcmp(name,"p")) { color = kColorProton; nProtons++;}
424 if (!strcmp(name,"p bar")) { color = kColorAntiProton; nAntiProtons++;}
425 if (!strcmp(name,"pi+")) { color = kColorPionPlus; nPionPlus++;}
426 if (!strcmp(name,"pi-")) { color = kColorPionMinus; nPionMinus++;}
427 if (!strcmp(name,"e+")) { color = kColorElectrons; nElectrons++;}
428 if (!strcmp(name,"e-")) { color = kColorElectrons; nElectrons++;}
429 if (!strcmp(name,"gamma")) { color = kColorGamma; nGammas++; lstyle = 3; }
430 if ( strstr(name,"K")) { color = kColorKaons; nKaons++;}
431 p->SetLineColor(color);
432 p->SetLineStyle(lstyle);
433 p->SetLineWidth(lwidth);
434 p->AppendPad();
435 ninvol++;
436 }
437
438 // event title
439 TPaveText *pt = new TPaveText(-0.94,0.85,-0.25,0.98,"br");
440 pt->AddText((char*)GetName());
441 pt->AddText((char*)GetTitle());
442 pt->SetFillColor(42);
443 pt->Draw();
444
445 // Annotate color codes
446 Int_t tcolor = 5;
447 if (gPad->GetFillColor() == 10) tcolor = 4;
448 TText *text = new TText(-0.95,-0.47,"Particles");
449 text->SetTextAlign(12);
450 text->SetTextSize(0.025);
451 text->SetTextColor(tcolor);
452 text->Draw();
453 text->SetTextColor(kColorGamma); text->DrawText(-0.95,-0.52,"(on screen)");
454 text->SetTextColor(kColorGamma); text->DrawText(-0.95,-0.57,"Gamma");
455 text->SetTextColor(kColorProton); text->DrawText(-0.95,-0.62,"Proton");
456 text->SetTextColor(kColorNeutron); text->DrawText(-0.95,-0.67,"Neutron");
457 text->SetTextColor(kColorAntiProton); text->DrawText(-0.95,-0.72,"AntiProton");
458 text->SetTextColor(kColorPionPlus); text->DrawText(-0.95,-0.77,"Pion +");
459 text->SetTextColor(kColorPionMinus); text->DrawText(-0.95,-0.82,"Pion -");
460 text->SetTextColor(kColorKaons); text->DrawText(-0.95,-0.87,"Kaons");
461 text->SetTextColor(kColorElectrons); text->DrawText(-0.95,-0.92,"Electrons,etc.");
462
463 text->SetTextColor(tcolor);
464 text->SetTextAlign(32);
465 char tcount[32];
466 snprintf(tcount,12,"%d",ntracks); text->DrawText(-0.55,-0.47,tcount);
467 snprintf(tcount,12,"%d",ninvol); text->DrawText(-0.55,-0.52,tcount);
468 snprintf(tcount,12,"%d",nGammas); text->DrawText(-0.55,-0.57,tcount);
469 snprintf(tcount,12,"%d",nProtons); text->DrawText(-0.55,-0.62,tcount);
470 snprintf(tcount,12,"%d",nNeutrons); text->DrawText(-0.55,-0.67,tcount);
471 snprintf(tcount,12,"%d",nAntiProtons); text->DrawText(-0.55,-0.72,tcount);
472 snprintf(tcount,12,"%d",nPionPlus); text->DrawText(-0.55,-0.77,tcount);
473 snprintf(tcount,12,"%d",nPionMinus); text->DrawText(-0.55,-0.82,tcount);
474 snprintf(tcount,12,"%d",nKaons); text->DrawText(-0.55,-0.87,tcount);
475 snprintf(tcount,12,"%d",nElectrons); text->DrawText(-0.55,-0.92,tcount);
476
477 text->SetTextAlign(12);
478 if (nPionPlus+nPionMinus) {
479 snprintf(tcount,31,"Protons/Pions= %4f",Float_t(nProtons)/Float_t(nPionPlus+nPionMinus));
480 } else {
481 strlcpy(tcount,"Protons/Pions= inf",31);
482 }
483 text->DrawText(-0.45,-0.92,tcount);
484
485 if (nPionPlus+nPionMinus) {
486 snprintf(tcount,31,"Kaons/Pions= %4f",Float_t(nKaons)/Float_t(nPionPlus+nPionMinus));
487 } else {
488 strlcpy(tcount,"Kaons/Pions= inf",31);
489 }
490 text->DrawText(0.30,-0.92,tcount);
491}
492
493
494////////////////////////////////////////////////////////////////////////////////
495/// Execute action corresponding to one event
496///
497
499{
500 if (gPad->GetView()) {
501 gPad->GetView()->ExecuteRotateView(event, px, py);
502 return;
503 }
504}
505
506////////////////////////////////////////////////////////////////////////////////
507/// Return the number of particles in the stack
508
510{
511 return fParticles->GetLast()+1;
512}
513
514////////////////////////////////////////////////////////////////////////////////
515/// Returns pointer to primary number i;
516///
517
519{
520 if (!fParticles) return 0;
522 if (i < 0 || i > n) return 0;
523 return (TParticle*)fParticles->UncheckedAt(i);
524}
525
526////////////////////////////////////////////////////////////////////////////////
527///
528/// Paint one event
529///
530
532{
533}
534
535////////////////////////////////////////////////////////////////////////////////
536///
537/// Set Pt threshold below which primaries are not drawn
538///
539
541{
542 fPtCut = ptcut;
543 Draw();
544 gPad->Update();
545}
546
547////////////////////////////////////////////////////////////////////////////////
548///
549/// Set lower and upper values of the view range
550///
551
553{
555}
556
557////////////////////////////////////////////////////////////////////////////////
558///
559/// Set lower and upper values of the view range
560///
561
563{
564 TView *view = gPad->GetView();
565 if (!view) return;
566 view->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
567
568 Draw();
569 gPad->Update();
570}
571
572////////////////////////////////////////////////////////////////////////////////
573///
574/// Set flag to display or not neutrons
575///
576
578{
579 fShowNeutrons = show;
580 Draw();
581 gPad->Update();
582}
#define HEPEVT
Definition Hepevt.h:32
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char text
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
#define gROOT
Definition TROOT.h:406
#define gPad
#define snprintf
Definition civetweb.c:1540
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
An array of clone (identical) objects.
void Clear(Option_t *option="") override
Clear the clones array.
static TDatabasePDG * Instance()
static function
The interface to various event generators.
Definition TGenerator.h:144
virtual TParticle * GetParticle(Int_t i) const
Returns pointer to primary number i;.
Int_t GetNumberOfParticles() const
Return the number of particles in the stack.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to objects in event.
virtual void SetViewRadius(Float_t rbox=1000)
Set lower and upper values of the view range.
~TGenerator() override
Event generator default destructor.
virtual void GenerateEvent()
must be implemented in concrete class (see eg TPythia6)
Bool_t fShowNeutrons
Pt cut. Do not show primaries below.
Definition TGenerator.h:148
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
It reads the /HEPEVT/ common block which has been filled by the GenerateEvent method.
virtual void SetPtCut(Float_t ptcut=0)
Set Pt threshold below which primaries are not drawn.
TObjArray * fParticles
display neutrons if true
Definition TGenerator.h:149
Float_t fPtCut
Definition TGenerator.h:147
void Browse(TBrowser *b) override
browse generator
void Draw(Option_t *option="") override
Insert one event in the pad list.
void Paint(Option_t *option="") override
Paint one event.
virtual void ShowNeutrons(Bool_t show=1)
Set flag to display or not neutrons.
virtual void SetViewRange(Float_t xmin=-10000, Float_t ymin=-10000, Float_t zmin=-10000, Float_t xmax=10000, Float_t ymax=10000, Float_t zmax=10000)
Set lower and upper values of the view range.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void Clear(Option_t *option="") override
Remove all objects from the array.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
Int_t GetLast() const override
Return index of last object in array.
void Add(TObject *obj) override
Definition TObjArray.h:68
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:202
Description of the static properties of a particle.
Description of the dynamic properties of a particle.
Definition TParticle.h:26
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
void Draw(Option_t *option="") override
Draw this pavetext with its current attributes.
Base class for several text objects.
Definition TText.h:22
See TView3D.
Definition TView.h:25
static TView * CreateView(Int_t system=1, const Double_t *rmin=nullptr, const Double_t *rmax=nullptr)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
virtual void GetRange(Float_t *min, Float_t *max)=0
virtual void SetRange(const Double_t *min, const Double_t *max)=0
TPaveText * pt
const Int_t n
Definition legend1.C:16
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:766
Definition rbox.py:1
th1 Draw()