FairRoot/PandaRoot
thailand2017/tut_fastsim.C
Go to the documentation of this file.
1 // *******
2 // Tutorial Macro for running fast simulation
3 // *******
4 
5 // The parameters are
6 // -------------------
7 // USAGE:\n";
8 // tut_fastsim.C+( [nevt], [prefix], [decfile], [mom], [res] )
9 // [nevt] : number of events; default = 1000
10 // [pref] : prefix for output
11 // [decfile] : decfile; 'DPM'/'FTF' uses DPM/FTF generator instead
12 // [mom] : pbar momentum
13 // [res] : resonance (ignored when running DPM); default = 'pbarpSystem0'
14 
15 
16 void tut_fastsim(Int_t nEvents = 1000, TString prefix = "signal", TString Decfile="pp_jpsi2pi_jpsi_mumu.dec", Float_t Mom=6.232, TString Resonance="pbarpSystem" )
17 {
18  // if Mom<0, interprete as -E_cm -> compute mom
19  double mp = 0.938272;
20  if (Mom<0)
21  {
22  double X = (Mom*Mom-2*mp*mp)/(2*mp);
23  Mom = sqrt(X*X-mp*mp);
24  }
25 
26  // Allow shortcut for resonance
27  if (Resonance=="pbp") Resonance = "pbarpSystem";
28  if (Resonance=="pbp0") Resonance = "pbarpSystem0";
29 
30  // Prevent generator from throwing a lot of warnings
31  TLorentzVector fIni(0,0,Mom,mp+sqrt(Mom*Mom+mp*mp));
32  TDatabasePDG::Instance()->AddParticle("pbarpSystem","pbarpSystem",fIni.M(),kFALSE,0.1,0, "",88888,0);
33  TDatabasePDG::Instance()->AddParticle("pbarpSystem0","pbarpSystem0",fIni.M(),kFALSE,0.1,0, "",88880,0);
34 
35  //----- Switches for Simulation Options ------------------------------
36  Bool_t enableSplitoff = true; // create e.-m. and hadronic split offs
37  Bool_t mergeNeutrals = true; // merge neutrals (for merged pi0s)
38  Bool_t electronBrems = true; // bremsstrahlung loss for electrons
39 
40  //----- Presist simulation output ------------------------------
41  Bool_t persist = true;
42 
43  //-----General settings-----------------------------------------------
44  TString BaseDir = gSystem->Getenv("VMCWORKDIR");
45  TString splitpars = BaseDir+"/fsim/splitpars.dat";
46  gRandom->SetSeed();
47 
48  //-----User Settings:-------------------------------------------------
49  TString OutputFile = prefix+"_fast.root";
50 
51  gDebug = 0;
52 
53  // choose your event generator
54  Bool_t UseEvtGenDirect = kTRUE;
55  Bool_t UseDpm = kFALSE;
56 
57  // use DPM generator; default: inelastic @ pbarmom = mom
58  if (Decfile.BeginsWith("DPM"))
59  {
60  UseEvtGenDirect = kFALSE;
61  UseDpm = kTRUE;
62  }
63 
64  // Start a stop watch
65  TStopwatch timer;
66  timer.Start();
67 
68  // Create the Simulation run manager
69  // --------------------------------
70  FairRunSim *fRun = new FairRunSim();
71  fRun->SetOutputFile(OutputFile.Data());
72  fRun->SetWriteRunInfoFile(kFALSE);
73 
74  FairLogger::GetLogger()->SetLogToFile(kFALSE);
75 
76  // -------------------------------
77  // Create and Set Event Generator
78  // -------------------------------
79  //FairFilteredPrimaryGenerator* primGen = new FairFilteredPrimaryGenerator();
80  FairPrimaryGenerator* primGen = new FairPrimaryGenerator();
81  fRun->SetGenerator(primGen);
82  fRun->SetName("TGeant3");
83 
84  if(UseDpm)
85  {
86  int mode = 0;
87  if (Decfile=="DPM1") mode = 1;
88  if (Decfile=="DPM2") mode = 2;
89 
90  PndDpmDirect *Dpm= new PndDpmDirect(Mom,mode); // 0 = inelastic, 1 = inelastic & elastic, 2 = elastic
91  Dpm->SetUnstable(111); // pi0
92  Dpm->SetUnstable(310); // K_S0
93  Dpm->SetUnstable(3122); // Lambda
94  Dpm->SetUnstable(-3122); // anti-Lambda
95  Dpm->SetUnstable(221); // eta
96  primGen->AddGenerator(Dpm);
97  }
98 
99  if(UseEvtGenDirect)
100  {
101  PndEvtGenDirect *EvtGen = new PndEvtGenDirect(Resonance, Decfile.Data(), Mom);
102  EvtGen->SetStoreTree(kTRUE);
103  primGen->AddGenerator(EvtGen);
104  }
105 
106  // ------------- switch off the transport of particles
107  primGen->DoTracking(kFALSE);
108 
109  //---------------------Create and Set the Field(s)----------
110  PndMultiField *fField= new PndMultiField("AUTO");
111  fRun->SetField(fField);
112 
113 
114  // ***********************************
115  // Setup the Fast Simulation Task
116  // ***********************************
117 
118  PndFastSim* fastSim = new PndFastSim(persist);
119 
120  // increasing verbosity increases the amount of console output (mainly for debugging)
121  fastSim->SetVerbosity(0);
122 
123  // enable the merging of neutrals if they have similar direction
124  //-----------------------------
125  fastSim->MergeNeutralClusters(mergeNeutrals);
126 
127  // enable bremsstahlung loss for electrons
128  //-----------------------------
129  fastSim->EnableElectronBremsstrahlung(electronBrems);
130 
131  //enable the producting of parametrized neutral (hadronic) split offs
132  // generate electro-magnetic / hadronic split offs in the EMC? switch off when running w/o EMC
133 
134  if (enableSplitoff) fastSim->EnableSplitoffs(splitpars.Data());
135 
136  fastSim->SetUseFlatCov(true);
137 
138  // -----------------------------------------------------------------------------------
139  // Tracking: Set up in parts of theta coverage. All modelled by PndFsmSimpleTracker.
140  // Mind: Numbers on resolution (pRes,thtRes,phiRes) and efficiency are guessed
141  // -----------------------------------------------------------------------------------
142 
143  // - (Full Panda Tracking: STT MVD GEM FTS)
144 
145  fastSim->AddDetector("ScSttAlone", "thtMin=145. thtMax=159.5 ptmin=0.1 pmin=0.0 pRes=0.04 thtRes=0.006 phiRes=0.007 efficiency=0.25");
146  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");
147  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");
148  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");
149 
150  // Fwd spectrometer enabled -> use Fwd tracking system
151  fastSim->AddDetector("ScFts", "thtMin=0. thtMax=5. ptmin=0.0 pmin=0.5 pRes=0.05 thtRes=0.002 phiRes=0.002 efficiency=0.80");
152 
153  // -----------------------------------------------------------------------------------
154  // Vertexing
155  // -----------------------------------------------------------------------------------
156 
157  // MVD and GEM are enabled -> better vertexing in central region
158  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
159  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
160 
161  // -----------------------------------------------------------------------------------
162  // EM Calorimeters w/ default parameters
163  // (don't have to be set, just to list the available parameters
164  // -----------------------------------------------------------------------------------
165 
166  fastSim->AddDetector("EmcFwCap", "thtMin=10.0 thtMax=22.0 Emin=0.01 dist=2.5");
167  fastSim->AddDetector("EmcBwCap", "thtMin=142.0 thtMax=160.0 Emin=0.01 dist=0.7");
168 
169  // EmcBarrel also allows to set phiMin and phiMax and can be added multiple times as EmcBarrel1, EmcBarrel2, etc.
170  // Should be made constistent with EmcPidBarrel below
171  fastSim->AddDetector("EmcBarrel","thtMin=22.0 thtMax=142.0 Emin=0.01 barrelRadius=0.5");
172 
173  // Fwd spectrometer enabled -> use Fwd EMC
174  fastSim->AddDetector("EmcFS", "thtMin=0.05 thtMax=10.0 aPar=0.013 bPar=0.0283 Emin=0.01 dist=8.2");
175 
176  // -----------------------------------------------------------------------------------
177  // PID
178  // -----------------------------------------------------------------------------------
179 
180  // PID detectors being always in: STT, MUO Barrel, EMC FwdCap, EMC BwdCap
181  //Note: A dEdX parametrization from 2008
182  fastSim->AddDetector("SttPid","thtMin=7.8 thtMax=159.5 ptmin=0.1 dEdxRes=0.15 efficiency=1.");
183  fastSim->AddDetector("ScMdtPidBarrel", "thtMin=10.0 thtMax=130.0 pmin=0.5 efficiency=0.95 misId=0.01");
184  fastSim->AddDetector("ScEmcPidFwCap", "thtMin=10.0 thtMax=22.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
185  fastSim->AddDetector("ScEmcPidBwCap", "thtMin=142.0 thtMax=160.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
186 
187  // MVD and GEM are enabled -> MVD PID available
188  //Note: A Bethe-Bloch-Landau-Gauss Prametrization from 2008
189  fastSim->AddDetector("MvdPid","thtMin=5. thtMax=133.6 ptmin=0.1 dEdxResMulti=1. efficiency=1.");
190 
191  // EMC Barrel enable -> EMC barrel PID available
192  fastSim->AddDetector("ScEmcPidBarrel", "thtMin=22.0 thtMax=142.0 ptmin=0.2 pmin=0.0 efficiency=1.0");
193 
194  // Barrel DIRC enabled
195  fastSim->AddDetector("DrcBarrel","thtMin=22.0 thtMax=140.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
196 
197  // Disc DIRC enabled
198  fastSim->AddDetector("DrcDisc","thtMin=5.0 thtMax=22.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
199 
200  // Fwd spectrometer enabled -> use RICH, FwdMUO and EMC FS
201  fastSim->AddDetector("ScEmcPidFS", "thtMin=0.5 thtMax=10.0 ptmin=0.0 pmin=0.5 efficiency=1.0");
202  fastSim->AddDetector("Rich","angleXMax=5.0 angleYMax=10.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
203  fastSim->AddDetector("ScMdtPidForward","thtMin=0.0 thtMax=10.0 pmin=0.5 efficiency=0.95 misId=0.01");
204 
205 
206  fRun->AddTask(fastSim);
207 
208  //------------------------- Initialize the RUN -----------------
209  fRun->Init();
210 
211  //------------------------- Run the Simulation -----------------
212  fRun->Run(nEvents);
213 
214 
215  //------------------------Print some info and exit----------------
216  timer.Stop();
217  Double_t rtime = timer.RealTime();
218  Double_t ctime = timer.CpuTime();
219  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
220 }
221 
PndMultiField * fField
Definition: sim_emc_apd.C:97
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
bool AddDetector(std::string name, std::string params="")
Definition: PndFastSim.cxx:313
static const double mp
Definition: mzparameters.h:11
void SetUseFlatCov(bool v=true)
Definition: PndFastSim.h:64
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
FairRunAna * fRun
Definition: hit_dirc.C:58
Int_t mode
Definition: autocutx.C:47
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
void tut_fastsim(Int_t nEvents=1000, TString prefix="signal", TString Decfile="pp_jpsi2pi_jpsi_mumu.dec", Float_t Mom=6.232, TString Resonance="pbarpSystem")
gDebug
Definition: sim_emc_apd.C:6
bool EnableSplitoffs(std::string fname="splitpars.dat")
Definition: PndFastSim.cxx:224
Double_t ctime
Definition: hit_dirc.C:114
double X
Definition: anaLmdDigi.C:68
void SetVerbosity(int vb)
Definition: PndFastSim.h:59
TLorentzVector fIni
Definition: full_core_ntp.C:29
Double_t rtime
Definition: hit_dirc.C:113