|
||||
File indexing completed on 2025-01-18 10:11:49
0001 // @(#)root/eg:$Id$ 0002 // Author: Ola Nordmann 21/09/95 0003 0004 /************************************************************************* 0005 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. * 0006 * All rights reserved. * 0007 * * 0008 * For the licensing terms see $ROOTSYS/LICENSE. * 0009 * For the list of contributors see $ROOTSYS/README/CREDITS. * 0010 *************************************************************************/ 0011 0012 0013 ////////////////////////////////////////////////////////////////////////// 0014 // // 0015 // TGenerator // 0016 // // 0017 // Is an base class, that defines the interface of ROOT to various // 0018 // event generators. Every event generator should inherit from // 0019 // TGenerator or its subclasses. // 0020 // // 0021 // Derived class can overload the member function GenerateEvent // 0022 // to do the actual event generation (e.g., call PYEVNT or similar). // 0023 // // 0024 // The derived class should overload the member function // 0025 // ImportParticles (both types) to read the internal storage of the // 0026 // generated event into either the internal TObjArray or the passed // 0027 // TClonesArray of TParticles. // 0028 // // 0029 // If the generator code stores event data in the /HEPEVT/ common block // 0030 // Then the default implementation of ImportParticles should suffice. // 0031 // The common block /HEPEVT/ is structed like // 0032 // // 0033 // /* C */ // 0034 // typedef struct { // 0035 // Int_t nevhep; // 0036 // Int_t nhep; // 0037 // Int_t isthep[4000]; // 0038 // Int_t idhep[4000]; // 0039 // Int_t jmohep[4000][2]; // 0040 // Int_t jdahep[4000][2]; // 0041 // Double_t phep[4000][5]; // 0042 // Double_t vhep[4000][4]; // 0043 // } HEPEVT_DEF; // 0044 // // 0045 // // 0046 // C Fortran // 0047 // COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(4000),IDHEP(4000), // 0048 // + JMOHEP(2,4000),JDAHEP(2,4000),PHEP(5,4000),VHEP(4,4000) // 0049 // INTEGER NEVHEP,NHEP,ISTHEP,IDHEP,JMOHEP,JDAHEP // 0050 // DOUBLE PRECISION PHEP,VHEP // 0051 // // 0052 // The generic member functions SetParameter and GetParameter can be // 0053 // overloaded to set and get parameters of the event generator. // 0054 // // 0055 // Note, if the derived class interfaces a (set of) Fortran common // 0056 // blocks (like TPythia, TVenus does), one better make the derived // 0057 // class a singleton. That is, something like // 0058 // // 0059 // class MyGenerator : public TGenerator // 0060 // { // 0061 // public: // 0062 // static MyGenerator* Instance() // 0063 // { // 0064 // if (!fgInstance) fgInstance = new MyGenerator; // 0065 // return fgInstance; // 0066 // } // 0067 // void GenerateEvent() { ... } // 0068 // void ImportParticles(TClonesArray* a, Option_t opt="") {...} // 0069 // Int_t ImportParticles(Option_t opt="") { ... } // 0070 // Int_t SetParameter(const char* name, Double_t val) { ... } // 0071 // Double_t GetParameter(const char* name) { ... } // 0072 // virtual ~MyGenerator() { ... } // 0073 // protected: // 0074 // MyGenerator() { ... } // 0075 // MyGenerator(const MyGenerator& o) { ... } // 0076 // MyGenerator& operator=(const MyGenerator& o) { ... } // 0077 // static MyGenerator* fgInstance; // 0078 // ClassDefOverride(MyGenerator,0); // 0079 // }; // 0080 // // 0081 // Having multiple objects accessing the same common blocks is not // 0082 // safe. // 0083 // // 0084 // concrete TGenerator classes can be loaded in scripts and subseqent- // 0085 // ly used in compiled code: // 0086 // // 0087 // // MyRun.h // 0088 // class MyRun : public TObject // 0089 // { // 0090 // public: // 0091 // static MyRun* Instance() { ... } // 0092 // void SetGenerator(TGenerator* g) { fGenerator = g; } // 0093 // void Run(Int_t n, Option_t* option="") // 0094 // { // 0095 // TFile* file = TFile::Open("file.root","RECREATE"); // 0096 // TTree* tree = new TTree("T","T"); // 0097 // TClonesArray* p = new TClonesArray("TParticles"); // 0098 // tree->Branch("particles", &p); // 0099 // for (Int_t event = 0; event < n; event++) { // 0100 // fGenerator->GenerateEvent(); // 0101 // fGenerator->ImportParticles(p,option); // 0102 // tree->Fill(); // 0103 // } // 0104 // file->Write(); // 0105 // file->Close(); // 0106 // } // 0107 // ... // 0108 // protected: // 0109 // TGenerator* fGenerator; // 0110 // ClassDefOverride(MyRun,0); // 0111 // }; // 0112 // // 0113 // // Config.C // 0114 // void Config() // 0115 // { // 0116 // MyRun* run = MyRun::Instance(); // 0117 // run->SetGenerator(MyGenerator::Instance()); // 0118 // } // 0119 // // 0120 // // main.cxx // 0121 // int // 0122 // main(int argc, char** argv) // 0123 // { // 0124 // TApplication app("", 0, 0); // 0125 // gSystem->ProcessLine(".x Config.C"); // 0126 // MyRun::Instance()->Run(10); // 0127 // return 0; // 0128 // } // 0129 // // 0130 // This is especially useful for example with TVirtualMC or similar. // 0131 // // 0132 ////////////////////////////////////////////////////////////////////////// 0133 0134 #ifndef ROOT_TGenerator 0135 #define ROOT_TGenerator 0136 0137 #include "TNamed.h" 0138 0139 class TBrowser; 0140 class TParticle; 0141 class TClonesArray; 0142 class TObjArray; 0143 0144 class TGenerator : public TNamed { 0145 0146 protected: 0147 Float_t fPtCut; //!Pt cut. Do not show primaries below 0148 Bool_t fShowNeutrons; //!display neutrons if true 0149 TObjArray *fParticles; //->static container of the primary particles 0150 0151 TGenerator(const TGenerator& tg) : 0152 TNamed(tg), fPtCut(tg.fPtCut), fShowNeutrons(tg.fShowNeutrons),fParticles(tg.fParticles) { } 0153 TGenerator& operator=(const TGenerator& tg) { 0154 if(this!=&tg) { 0155 TNamed::operator=(tg); fPtCut=tg.fPtCut; fShowNeutrons=tg.fShowNeutrons; 0156 fParticles=tg.fParticles; 0157 } 0158 return *this; 0159 } 0160 0161 public: 0162 0163 TGenerator(): fPtCut(0), fShowNeutrons(kTRUE), fParticles(0) { } //Used by Dictionary 0164 TGenerator(const char *name, const char *title="Generator class"); 0165 ~TGenerator() override; 0166 void Browse(TBrowser *b) override; 0167 Int_t DistancetoPrimitive(Int_t px, Int_t py) override; 0168 void Draw(Option_t *option="") override; 0169 void ExecuteEvent(Int_t event, Int_t px, Int_t py) override; 0170 virtual void GenerateEvent(); 0171 virtual Double_t GetParameter(const char* /*name*/) const { return 0.; } 0172 virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option=""); 0173 virtual TObjArray *ImportParticles(Option_t *option=""); 0174 virtual TParticle *GetParticle(Int_t i) const; 0175 Int_t GetNumberOfParticles() const; 0176 virtual TObjArray *GetListOfParticles() const {return fParticles;} 0177 virtual TObjArray *GetPrimaries(Option_t *option="") {return ImportParticles(option);} 0178 Float_t GetPtCut() const {return fPtCut;} 0179 void Paint(Option_t *option="") override; 0180 virtual void SetParameter(const char* /*name*/,Double_t /*val*/){} 0181 virtual void SetPtCut(Float_t ptcut=0); // *MENU* 0182 virtual void SetViewRadius(Float_t rbox = 1000); // *MENU* 0183 virtual void SetViewRange(Float_t xmin=-10000,Float_t ymin=-10000,Float_t zmin=-10000 0184 ,Float_t xmax=10000,Float_t ymax=10000,Float_t zmax=10000); // *MENU* 0185 virtual void ShowNeutrons(Bool_t show=1); // *MENU* 0186 0187 ClassDefOverride(TGenerator,1); //Event generator interface abstract baseclass 0188 }; 0189 0190 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |