Logo ROOT   6.14/05
Reference Guide
TGenerator.h
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 
13 //////////////////////////////////////////////////////////////////////////
14 // //
15 // TGenerator //
16 // //
17 // Is an base class, that defines the interface of ROOT to various //
18 // event generators. Every event generator should inherit from //
19 // TGenerator or its subclasses. //
20 // //
21 // Derived class can overload the member function GenerateEvent //
22 // to do the actual event generation (e.g., call PYEVNT or similar). //
23 // //
24 // The derived class should overload the member function //
25 // ImportParticles (both types) to read the internal storage of the //
26 // generated event into either the internal TObjArray or the passed //
27 // TClonesArray of TParticles. //
28 // //
29 // If the generator code stores event data in the /HEPEVT/ common block //
30 // Then the default implementation of ImportParticles should suffice. //
31 // The common block /HEPEVT/ is structed like //
32 // //
33 // /* C */ //
34 // typedef struct { //
35 // Int_t nevhep; //
36 // Int_t nhep; //
37 // Int_t isthep[4000]; //
38 // Int_t idhep[4000]; //
39 // Int_t jmohep[4000][2]; //
40 // Int_t jdahep[4000][2]; //
41 // Double_t phep[4000][5]; //
42 // Double_t vhep[4000][4]; //
43 // } HEPEVT_DEF; //
44 // //
45 // //
46 // C Fortran //
47 // COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(4000),IDHEP(4000), //
48 // + JMOHEP(2,4000),JDAHEP(2,4000),PHEP(5,4000),VHEP(4,4000) //
49 // INTEGER NEVHEP,NHEP,ISTHEP,IDHEP,JMOHEP,JDAHEP //
50 // DOUBLE PRECISION PHEP,VHEP //
51 // //
52 // The generic member functions SetParameter and GetParameter can be //
53 // overloaded to set and get parameters of the event generator. //
54 // //
55 // Note, if the derived class interfaces a (set of) Fortran common //
56 // blocks (like TPythia, TVenus does), one better make the derived //
57 // class a singleton. That is, something like //
58 // //
59 // class MyGenerator : public TGenerator //
60 // { //
61 // public: //
62 // static MyGenerator* Instance() //
63 // { //
64 // if (!fgInstance) fgInstance = new MyGenerator; //
65 // return fgInstance; //
66 // } //
67 // void GenerateEvent() { ... } //
68 // void ImportParticles(TClonesArray* a, Option_t opt="") {...} //
69 // Int_t ImportParticles(Option_t opt="") { ... } //
70 // Int_t SetParameter(const char* name, Double_t val) { ... } //
71 // Double_t GetParameter(const char* name) { ... } //
72 // virtual ~MyGenerator() { ... } //
73 // protected: //
74 // MyGenerator() { ... } //
75 // MyGenerator(const MyGenerator& o) { ... } //
76 // MyGenerator& operator=(const MyGenerator& o) { ... } //
77 // static MyGenerator* fgInstance; //
78 // ClassDef(MyGenerator,0); //
79 // }; //
80 // //
81 // Having multiple objects accessing the same common blocks is not //
82 // safe. //
83 // //
84 // concrete TGenerator classes can be loaded in scripts and subseqent- //
85 // ly used in compiled code: //
86 // //
87 // // MyRun.h //
88 // class MyRun : public TObject //
89 // { //
90 // public: //
91 // static MyRun* Instance() { ... } //
92 // void SetGenerator(TGenerator* g) { fGenerator = g; } //
93 // void Run(Int_t n, Option_t* option="") //
94 // { //
95 // TFile* file = TFile::Open("file.root","RECREATE"); //
96 // TTree* tree = new TTree("T","T"); //
97 // TClonesArray* p = new TClonesArray("TParticles"); //
98 // tree->Branch("particles", &p); //
99 // for (Int_t event = 0; event < n; event++) { //
100 // fGenerator->GenerateEvent(); //
101 // fGenerator->ImportParticles(p,option); //
102 // tree->Fill(); //
103 // } //
104 // file->Write(); //
105 // file->Close(); //
106 // } //
107 // ... //
108 // protected: //
109 // TGenerator* fGenerator; //
110 // ClassDef(MyRun,0); //
111 // }; //
112 // //
113 // // Config.C //
114 // void Config() //
115 // { //
116 // MyRun* run = MyRun::Instance(); //
117 // run->SetGenerator(MyGenerator::Instance()); //
118 // } //
119 // //
120 // // main.cxx //
121 // int //
122 // main(int argc, char** argv) //
123 // { //
124 // TApplication app("", 0, 0); //
125 // gSystem->ProcessLine(".x Config.C"); //
126 // MyRun::Instance()->Run(10); //
127 // return 0; //
128 // } //
129 // //
130 // This is especially useful for example with TVirtualMC or similar. //
131 // //
132 //////////////////////////////////////////////////////////////////////////
133 
134 #ifndef ROOT_TGenerator
135 #define ROOT_TGenerator
136 
137 #include "TNamed.h"
138 
139 class TBrowser;
140 class TParticle;
141 class TClonesArray;
142 class TObjArray;
143 
144 class TGenerator : public TNamed {
145 
146 protected:
147  Float_t fPtCut; //!Pt cut. Do not show primaries below
148  Bool_t fShowNeutrons; //!display neutrons if true
149  TObjArray *fParticles; //->static container of the primary particles
150 
151  TGenerator(const TGenerator& tg) :
152  TNamed(tg), fPtCut(tg.fPtCut), fShowNeutrons(tg.fShowNeutrons),fParticles(tg.fParticles) { }
154  if(this!=&tg) {
155  TNamed::operator=(tg); fPtCut=tg.fPtCut; fShowNeutrons=tg.fShowNeutrons;
156  fParticles=tg.fParticles;
157  }
158  return *this;
159  }
160 
161 public:
162 
163  TGenerator(): fPtCut(0), fShowNeutrons(kTRUE), fParticles(0) { } //Used by Dictionary
164  TGenerator(const char *name, const char *title="Generator class");
165  virtual ~TGenerator();
166  virtual void Browse(TBrowser *b);
167  virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
168  virtual void Draw(Option_t *option="");
169  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
170  virtual void GenerateEvent();
171  virtual Double_t GetParameter(const char* /*name*/) const { return 0.; }
172  virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="");
173  virtual TObjArray *ImportParticles(Option_t *option="");
174  virtual TParticle *GetParticle(Int_t i) const;
175  Int_t GetNumberOfParticles() const;
176  virtual TObjArray *GetListOfParticles() const {return fParticles;}
177  virtual TObjArray *GetPrimaries(Option_t *option="") {return ImportParticles(option);}
178  Float_t GetPtCut() const {return fPtCut;}
179  virtual void Paint(Option_t *option="");
180  virtual void SetParameter(const char* /*name*/,Double_t /*val*/){}
181  virtual void SetPtCut(Float_t ptcut=0); // *MENU*
182  virtual void SetViewRadius(Float_t rbox = 1000); // *MENU*
183  virtual void SetViewRange(Float_t xmin=-10000,Float_t ymin=-10000,Float_t zmin=-10000
184  ,Float_t xmax=10000,Float_t ymax=10000,Float_t zmax=10000); // *MENU*
185  virtual void ShowNeutrons(Bool_t show=1); // *MENU*
186 
187  ClassDef(TGenerator,1) //Event generator interface abstract baseclass
188 };
189 
190 #endif
An array of TObjects.
Definition: TObjArray.h:37
float xmin
Definition: THbookFile.cxx:93
Float_t fPtCut
Definition: TGenerator.h:147
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Bool_t fShowNeutrons
Pt cut. Do not show primaries below.
Definition: TGenerator.h:148
float ymin
Definition: THbookFile.cxx:93
The interface to various event generators.
Definition: TGenerator.h:144
Int_t GetNumberOfParticles() const
Return the number of particles in the stack.
Definition: TGenerator.cxx:506
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
TGenerator & operator=(const TGenerator &tg)
Definition: TGenerator.h:153
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TGenerator.cxx:495
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
It reads the /HEPEVT/ common block which has been filled by the GenerateEvent method.
Definition: TGenerator.cxx:270
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
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.
Definition: TGenerator.cxx:559
virtual Double_t GetParameter(const char *) const
Definition: TGenerator.h:171
Float_t GetPtCut() const
Definition: TGenerator.h:178
#define ClassDef(name, id)
Definition: Rtypes.h:320
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual void ShowNeutrons(Bool_t show=1)
Set flag to display or not neutrons.
Definition: TGenerator.cxx:574
virtual void GenerateEvent()
must be implemented in concrete class (see eg TPythia6)
Definition: TGenerator.cxx:189
virtual TObjArray * GetPrimaries(Option_t *option="")
Definition: TGenerator.h:177
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
virtual void SetParameter(const char *, Double_t)
Definition: TGenerator.h:180
float ymax
Definition: THbookFile.cxx:93
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition: TNamed.cxx:51
virtual void Draw(Option_t *option="")
Insert one event in the pad list.
Definition: TGenerator.cxx:348
virtual void Paint(Option_t *option="")
Paint one event.
Definition: TGenerator.cxx:528
TGenerator(const TGenerator &tg)
Definition: TGenerator.h:151
float xmax
Definition: THbookFile.cxx:93
virtual void SetViewRadius(Float_t rbox=1000)
Set lower and upper values of the view range.
Definition: TGenerator.cxx:549
double Double_t
Definition: RtypesCore.h:55
TObjArray * fParticles
display neutrons if true
Definition: TGenerator.h:149
virtual void SetPtCut(Float_t ptcut=0)
Set Pt threshold below which primaries are not drawn.
Definition: TGenerator.cxx:537
An array of clone (identical) objects.
Definition: TClonesArray.h:32
virtual void Browse(TBrowser *b)
browse generator
Definition: TGenerator.cxx:324
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
virtual ~TGenerator()
Event generator default destructor.
Definition: TGenerator.cxx:176
virtual TParticle * GetParticle(Int_t i) const
Returns pointer to primary number i;.
Definition: TGenerator.cxx:515
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual TObjArray * GetListOfParticles() const
Definition: TGenerator.h:176
char name[80]
Definition: TGX11.cxx:109
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to objects in event.
Definition: TGenerator.cxx:334