FairRoot/PandaRoot
macro/production/scrut-2014/simfast.C
Go to the documentation of this file.
1 // *******
2 // Macro for running fast simulation
3 // *******
4 
5 // The parameters are
6 // -------------------
7 // Prefix : Prefix string for the output file
8 // Decfile : 1) name of the EvtGen decay file; 2) 'DPM' for using DpmGen; 3) 'BOX' for using Box Generator
9 // Mom : 1) EvtGen and DpmGen: pbar momentum in GeV/c; 2) BoxGen, Mom > 0: generate 0.1 < p < Mom; 3) BoxGen, Mom < 0: generate only p = -Mom
10 // nEvents : Number of events to be generated
11 // Resonance : Initial resonance name for EvtGen
12 // Pdgcode : Only BoxGen: pdgcode of particle to be generated
13 
14 int simfast(TString Prefix="test", TString Decfile="psi2s_Jpsi2pi_Jpsi_mumu.dec", Float_t Mom=6.231552, Int_t nEvents = 10000, TString Resonance="pbarpSystem", int Pdgcode = 11 )
15 {
16  TString BaseDir = gSystem->Getenv("VMCWORKDIR");
17  TString splitpars = BaseDir+"/fsim/splitpars.dat";
18  gRandom->SetSeed();
19 
20  //-----User Settings:-----------------------------------------------
21  TString OutputFile = Prefix+"_fast.root";
22  gDebug = 0;
23 
24  // generate electro-magnetic / hadronic split offs in the EMC? switch off when running w/o EMC
25  Bool_t enableSplitoff = kFALSE;
26 
27  // choose your event generator
28  Bool_t UseEvtGenDirect = kTRUE;
29  Bool_t UseDpm = kFALSE;
30  Bool_t UseBoxGenerator = kFALSE;
31 
32  // use DPM generator; default: inelastic @ pbarmom = mom
33  if (Decfile=="DPM")
34  {
35  UseEvtGenDirect = kFALSE;
36  UseDpm = kTRUE;
37  }
38 
39  // use BOX generator; default: single mu-, 0<tht<180, 0<phi<360, 0.1<p<mom
40  if (Decfile=="BOX")
41  {
42  UseEvtGenDirect = kFALSE;
43  UseBoxGenerator = kTRUE;
44  }
45 
46  Double_t MomMin = 0.1; // minimum momentum for box generator
47  Double_t MomMax = Mom; // maximum " "
48 
49  // Start a stop watch
50  TStopwatch timer;
51  timer.Start();
52 
53  // Create the Simulation run manager
54  // --------------------------------
55  FairRunSim *fRun = new FairRunSim();
56  fRun->SetOutputFile(OutputFile.Data());
57  fRun->SetGenerateRunInfo(kFALSE);
58 
59  FairLogger::GetLogger()->SetLogToFile(kFALSE);
60 
61 
62  // Create and Set Event Generator
63  // -------------------------------
64  FairPrimaryGenerator* primGen = new FairPrimaryGenerator();
65  fRun->SetGenerator(primGen);
66  fRun->SetName("TGeant3");
67 
68  if(UseBoxGenerator){ // Box Generator
69  FairBoxGenerator* boxGen = new FairBoxGenerator(Pdgcode, 1); // 211 = pion; 1 = multipl.
70  boxGen->SetPRange(MomMin,MomMax); // GeV/c
71  boxGen->SetPhiRange(0., 360.); // Azimuth angle range [degree]
72  boxGen->SetThetaRange(0., 180.); // Polar angle in lab system range [degree]
73  boxGen->SetXYZ(0., 0., 0.); //cm
74  primGen->AddGenerator(boxGen);
75  }
76  if(UseDpm){
77  PndDpmDirect *Dpm= new PndDpmDirect(Mom,0); // 0 = inelastic, 1 = inelastic & elastic, 2 = elastic
78  Dpm->SetUnstable(111); // pi0
79  Dpm->SetUnstable(310); // K_S0
80  Dpm->SetUnstable(3122); // Lambda
81  Dpm->SetUnstable(-3122); // anti-Lambda
82  Dpm->SetUnstable(221); // eta
83  primGen->AddGenerator(Dpm);
84  }
85  if(UseEvtGenDirect){
86  PndEvtGenDirect *EvtGen = new PndEvtGenDirect(Resonance, Decfile.Data(), Mom);
87  EvtGen->SetStoreTree(kTRUE);
88  primGen->AddGenerator(EvtGen);
89  }
90 
91  // ------------- switch off the transport of particles
92  primGen->DoTracking(kFALSE);
93 
94  //---------------------Create and Set the Field(s)----------
95  PndMultiField *fField= new PndMultiField("FULL");
96  fRun->SetField(fField);
97 
98  //Setup the Fast Simulation Task
99  //-----------------------------
100  PndFastSim* fastSim = new PndFastSim();
101 
102  //increasing verbosity increases the amount of console output (mainly for debugging)
103  fastSim->SetVerbosity(0);
104 
105  // enable the merging of neutrals if they have similar direction
106  fastSim->MergeNeutralClusters(false);
107 
108  // enable bremsstahlung loss for electrons
109  fastSim->EnableElectronBremsstrahlung(false);
110 
111  //enable the producting of parametrized neutral (hadronic) split offs
112  if (enableSplitoff)
113  fastSim->EnableSplitoffs(splitpars.Data());
114 
115  fastSim->SetUseFlatCov(true);
116  // -----------------------------------------------------------------------------------
117  //Tracking: Set up in parts of theta coverage. All modelled by PndFsmSimpleTracker.
118  // Mind: Numbers on resolution (pRes,thtRes,phiRes) and efficiency are guessed
119  // -----------------------------------------------------------------------------------
120  // - (Full Panda Tracking: STT MVD GEM FTS)
121  fastSim->AddDetector("ScSttAlone", "thtMin=145. thtMax=159.5 ptmin=0.1 pmin=0.0 pRes=0.04 thtRes=0.001 phiRes=0.001 efficiency=0.25");
122  fastSim->AddDetector("ScSttMvd", "thtMin=20.9 thtMax=145. ptmin=0.1 pmin=0.0 pRes=0.02 thtRes=0.001 phiRes=0.001 efficiency=0.85");
123  fastSim->AddDetector("ScSttMvdGem", "thtMin=7.8 thtMax=20.9 ptmin=0.1 pmin=0.0 pRes=0.02 thtRes=0.001 phiRes=0.001 efficiency=0.85");
124  fastSim->AddDetector("ScMvdGem", "thtMin=5. thtMax=7.8 ptmin=0.1 pmin=0.0 pRes=0.03 thtRes=0.001 phiRes=0.001 efficiency=0.60");
125  fastSim->AddDetector("ScFts", "thtMin=0.0 thtMax=5. ptmin=0.0 pmin=0.5 pRes=0.05 thtRes=0.002 phiRes=0.002 efficiency=0.80");
126  // - other options:
127  //fastSim->AddDetector("ScMvdGemFts", "thtMin=5. thtMax=7.8 ptmin=0.1 pmin=0.0 pRes=0.03 thtRes=0.001 phiRes=0.001 efficiency=0.80");
128  //fastSim->AddDetector("ScMvdFts", "thtMin=5. thtMax=10 ptmin=0.0 pmin=0.0 pRes=0.05 thtRes=0.002 phiRes=0.002 efficiency=0.80");
129  //fastSim->AddDetector("ScGemFts", "thtMin=3. thtMax=10 ptmin=0.0 pmin=0.0 pRes=0.05 thtRes=0.002 phiRes=0.002 efficiency=0.80");
130  // - STT alone:
131  //fastSim->AddDetector("ScSttAlone", "thtMin=133.6 thtMax=159.5 ptmin=0.1 pmin=0.0 pRes=0.030 thtRes=0.06 phiRes=0.1 efficiency=0.25");
132  //fastSim->AddDetector("ScSttAlone2", "thtMin=20.9 thtMax=133.6 ptmin=0.1 pmin=0.0 pRes=0.026 thtRes=0.06 phiRes=0.1 efficiency=0.95");
133  //fastSim->AddDetector("ScSttAlone3", "thtMin=7.8 thtMax=20.9 ptmin=0.1 pmin=0.0 pRes=0.026 thtRes=0.06 phiRes=0.1 efficiency=0.25");
134 
135  // -----------------------------------------------------------------------------------
136  // Vertexing
137  // -----------------------------------------------------------------------------------
138  // - A
139  fastSim->AddDetector("ScVtxMvd", "thtMin=5. thtMax=145. ptmin=0.1 vtxRes=0.005 efficiency=1."); // efficiency=1: all tracks found in trackers will get a vertex information
140  fastSim->AddDetector("ScVtxNoMvd", "thtMin=0. thtMax=5. ptmin=0.0 vtxRes=0.05 efficiency=1."); // efficiency=1: all tracks found in trackers will get a vertex information
141  // - B
142  //fastSim->AddDetector("ScVtxNoMvd", "thtMin=0. thtMax=160. ptmin=0.1 vtxRes=0.1 efficiency=1."); // efficiency=1: all tracks found in trackers will get a vertex information
143 
144  // -----------------------------------------------------------------------------------
145  // EM Calorimeters w/ default parameters
146  // (don't have to be set, just to list the available parameters
147  // EmcBarrel also allows to set phiMin and phiMax and can be added multiple times as EmcBarrel1, EmcBarrel2, etc.
148  // Should be made constistent with EmcPidBarrel below
149  // -----------------------------------------------------------------------------------
150  fastSim->AddDetector("EmcBarrel","thtMin=22.0 thtMax=142.0 Emin=0.01 barrelRadius=0.5");
151  fastSim->AddDetector("EmcFwCap", "thtMin=10.0 thtMax=22.0 Emin=0.01 dist=2.5");
152  fastSim->AddDetector("EmcBwCap", "thtMin=142.0 thtMax=160.0 Emin=0.01 dist=0.7");
153  fastSim->AddDetector("EmcFS", "thtMin=0.05 thtMax=10.0 aPar=0.02 bPar=0.0274 Emin=0.01 dist=8.0");
154 
155  // -----------------------------------------------------------------------------------
156  // PID
157  //
158  // -----------------------------------------------------------------------------------
159  fastSim->AddDetector("MvdPid","thtMin=5. thtMax=133.6 ptmin=0.1 dEdxResMulti=1. efficiency=1."); //Note: A Bethe-Bloch-Landau-Gauss Prametrization from 2008
160  fastSim->AddDetector("SttPid","thtMin=7.8 thtMax=159.5 ptmin=0.1 dEdxRes=1. efficiency=1."); //Note: A dEdX parametrization from 2008
161 
162  fastSim->AddDetector("DrcBarrel","thtMin=22.0 thtMax=140.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
163  fastSim->AddDetector("DrcDisc","thtMin=5.0 thtMax=22.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
164  fastSim->AddDetector("Rich","angleXMax=5.0 angleYMax=10.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
165 
166  fastSim->AddDetector("ScMdtPidBarrel", "thtMin=10.0 thtMax=130.0 pmin=0.5 efficiency=0.95 misId=0.01");
167  fastSim->AddDetector("ScMdtPidForward","thtMin=0.0 thtMax=10.0 pmin=0.5 efficiency=0.95 misId=0.01");
168 
169  fastSim->AddDetector("ScEmcPidFS", "thtMin=0.5 thtMax=10.0 ptmin=0.0 pmin=0.5 efficiency=1.0");
170  fastSim->AddDetector("ScEmcPidFwCap", "thtMin=10.0 thtMax=22.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
171  fastSim->AddDetector("ScEmcPidBarrel", "thtMin=22.0 thtMax=142.0 ptmin=0.2 pmin=0.0 efficiency=1.0");
172  fastSim->AddDetector("ScEmcPidBwCap", "thtMin=142.0 thtMax=160.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
173 
174  fRun->AddTask(fastSim);
175  //------------------------- Initialize the RUN -----------------
176  fRun->Init();
177  //------------------------- Run the Simulation -----------------
178  fRun->Run(nEvents);
179  //------------------------Print some info and exit----------------
180  timer.Stop();
181  Double_t rtime = timer.RealTime();
182  Double_t ctime = timer.CpuTime();
183  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
184  return 0;
185 }
186 
PndMultiField * fField
Definition: sim_emc_apd.C:97
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
bool AddDetector(std::string name, std::string params="")
Definition: PndFastSim.cxx:313
void SetUseFlatCov(bool v=true)
Definition: PndFastSim.h:64
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
FairRunAna * fRun
Definition: hit_dirc.C:58
void SetUnstable(int pdg)
bool MergeNeutralClusters(bool merge=true, double par=0.389)
Definition: PndFastSim.h:61
void SetStoreTree(Bool_t store=true)
Double_t
void EnableElectronBremsstrahlung(bool brems=true)
Definition: PndFastSim.h:62
Int_t nEvents
Definition: hit_dirc.C:11
TStopwatch timer
Definition: hit_dirc.C:51
gDebug
Definition: sim_emc_apd.C:6
int simfast(TString Prefix="test", TString Decfile="psi2s_Jpsi2pi_Jpsi_mumu.dec", Float_t Mom=6.231552, Int_t nEvents=10000, TString Resonance="pbarpSystem", int Pdgcode=11)
bool EnableSplitoffs(std::string fname="splitpars.dat")
Definition: PndFastSim.cxx:224
Double_t ctime
Definition: hit_dirc.C:114
FairBoxGenerator * boxGen
Definition: sim_emc_apd.C:85
void SetVerbosity(int vb)
Definition: PndFastSim.h:59
Double_t rtime
Definition: hit_dirc.C:113