FairRoot/PandaRoot
macro/production/dayone-2017/fastsim/quickfsimana.C
Go to the documentation of this file.
1 // **********************************************************************************************
2 // Macro for running fast simulation + Software Trigger + simple analysis (PndSimpleCombinerTask)
3 // **********************************************************************************************
4 
5 // The parameters are
6 // -------------------
7 // quickfsimana.C( <pref>, <decfile>, <mom>, <decay>, [nevt], [res], [parms], [runST], [runnum], [mode] )
8 // <pref> : output file names prefix
9 // <decfile> : EvtGen decfile; DPM/FTF/BOX uses DPM/FTF generator (inelastic mode) or box generator instead
10 // <mom> : EvtGen, DPM, FTF: pbar momentum (negative values are interpreted as -E_cm); BOX generator: maximum particle momentum
11 // [decay] : the decay pattern to be reconstructed, e.g. 'phi -> K+ K-; D_s+ -> phi pi+'; emtpy string: only fast sim will be run
12 // [nevt] : number of events; default = 1000
13 // [res] : initial resonance or particle type for BOX generator (ignored when running DPM); default = 'pbarpSystem0'
14 // [parms] : parameters for the analysis, e.g. 'mwin=0.4:mwin(phi)=0.1:emin=0.1:pmin=0.1:qamc'; 'qapart' runs PndParticleQATask
15 // [runST] : if 'true' runs Software Trigger (default: false)
16 // [runnum] : integer run number (default: 0)
17 // [mode] : arbitrary mode number (default: 0)
18 // -------------------
19 
20 void getRange(TString par, double &min, double &max);
21 TString getInitialResonance(TString &fEvtGenFile);
22 
23 void quickfsimana(TString Prefix="", TString Decfile="", Float_t Mom=0., TString anadecay="", Int_t nEvents = 1000, TString anaparms="", bool runST=false, int run=0 , int runmode=0)
24 {
25  if (Prefix=="" || Decfile=="" || Mom==0. )
26  {
27  cout << "USAGE:\n";
28  cout << "quickfsimana.C( <pref>, <decfile>, <mom>, <decay>, [nevt], [parms], [runST], [runnum], [mode] )\n\n";
29  cout << " <pref> : output file names prefix\n";
30  cout << " <decfile> : EvtGen decfile 'xxx.dec' or 'xxx.dec:iniRes'; DPM/FTF/BOX uses DPM/FTF generator or BOX generator instead\n";
31  cout << " DPM settings: DPM = inelastic only, DPM1 = inel. + elastic, DPM2 = elastic only\n";
32  cout << " FTF settings: FTF = inelastic only, FTF1 = inel. + elastic\n";
33  cout << " BOX settings: optional ranges 'p/tht/phi[min,max]' separated with colon; single number = fixed value; example: 'BOX:p[1,5]:tht[45]:phi[90,210]'\n";
34  cout << " <mom> : EvtGen, DPM, FTF: pbar momentum (negative values are interpreted as -E_cm); BOX generator w/o special settings: maximum particle momentum\n";
35  cout << " [decay] : the decay pattern to be reconstructed, e.g. 'phi -> K+ K-; D_s+ -> phi pi+'; '': only fast sim w/o reco will be run\n";
36  cout << " [nevt] : number of events; default = 1000\n";
37  //cout << " [res] : initial resonance or particle type for BOX generator (ignored when running DPM); default = 'pbarpSystem0'\n";
38  cout << " [parms] : parameters for the analysis, e.g. 'mwin=0.4:mwin(phi)=0.1:emin=0.1:pmin=0.1:qamc'; 'qapart' runs PndParticleQATask: 'persist' saves PndPidCandidates\n";
39  cout << " [runST] : if 'true' runs Software Trigger (default: false)\n";
40  cout << " [runnum] : integer run number (default: 0)\n";
41  cout << " [mode] : arbitrary mode number (default: 0)\n\n";
42  cout << "Example 1 - Do reco for EvtGen events : root -l -b -q 'quickfsimana.C(\"jpsi2pi\", \"pp_jpsi2pi_jpsi_mumu.dec\", 6.23, \"J/psi -> mu+ mu-; pbp -> J/psi pi+ pi-\", 1000, \"fit4c:mwin=0.6\")'\n";
43  cout << "Example 2 - Particle QA for BOX gen : root -l -b -q 'quickfsimana.C(\"single_kplus\", \"BOX:type(321,1)\", 8.0, \"\", 1000, \"qapart\")'\n";
44  cout << "Example 3 - Run fast sim only for DPM : root -l -b -q 'quickfsimana.C(\"bkg\", \"DPM\", 6.23, \"\", 1000)'\n\n";
45  return;
46  }
47 
48  // persist fast sim output?
49  bool persist = (anadecay == "" && anaparms == "") || anaparms.Contains("persist");
50 
51  // do some reconstruction ?
52  bool doreco = (anadecay != "" || anaparms.Contains("nevt"));
53 
54  // do particle QA?
55  bool partQA = (anaparms.Contains("qapart"));
56  bool mc = !(anaparms.Contains("!mc")) && !(anaparms.Contains("qamc"));
57  bool neut = !(anaparms.Contains("!neut"));
58  bool chrg = !(anaparms.Contains("!chrg"));
59 
60  // for submission to queue all blanks in decay string were replaced by '§'; now we replace again the other way around
61  anadecay.ReplaceAll("#"," ");
62 
63  // if Mom<0, interprete as -E_cm
64  double mp = 0.938272;
65  if (Mom<0)
66  {
67  double X = (Mom*Mom-2*mp*mp)/(2*mp);
68  Mom = sqrt(X*X-mp*mp);
69  }
70 
71  // Allow shortcut for resonance pbarpSystem
72  anadecay.ReplaceAll("pbp", "pbarpSystem");
73 
74  // Prevent generator from throwing a lot of warnings
75  TLorentzVector fIni(0,0,Mom,mp+sqrt(Mom*Mom+mp*mp));
76  TDatabasePDG *pdg=TDatabasePDG::Instance();
77  pdg->AddParticle("pbarpSystem","pbarpSystem",fIni.M(),kFALSE,0.1,0, "",88888,0);
78  pdg->AddParticle("pbarpSystem0","pbarpSystem0",fIni.M(),kFALSE,0.1,0, "",88880,0);
79  pdg->AddParticle("eta_c1","eta_c1",4.3,kFALSE,0.02,0,"",999441);
80 
81  //-----Evaluate Detector Setup ---------------------------------------
82  bool SwMvdGem = true; // Enable MVD and GEM for central tracking in addition to STT
83  bool SwEmcBar = true; // Enable EMC barrel for calorimetry (neutral detection and PID component)
84  bool SwDrc = true; // Enable Barrel DIRC for PID
85  bool SwDsc = true; // Enable Disc DIRC for PID
86  bool SwFwdSpec = true; // Enable complete Forward Spectrometer (= Fwd Spec. EMC, Fwd Tracking, RICH, Fwd MUO)
87 
88  //----- Switches for Simulation Options ------------------------------
89  Bool_t enableSplitoff = true; // create e.-m. and hadronic split offs
90  Bool_t mergeNeutrals = true; // merge neutrals (for merged pi0s)
91  Bool_t electronBrems = true; // bremsstrahlung loss for electrons
92  Bool_t useEventFilter = false; // enable Fast Sim event filter. *** Needs configuration (see below) ***
93  Bool_t usePndEventFilter = false; // enable Panda event filter. *** Needs configuration (see below) ***
94 
95  //-----General settings-----------------------------------------------
96  TString BaseDir = gSystem->Getenv("VMCWORKDIR");
97  TString splitpars = BaseDir+"/fsim/splitpars.dat";
98  gRandom->SetSeed();
99 
100  //-----User Settings:-------------------------------------------------
101  TString OutputFile = TString::Format("%s_%d_ana.root",Prefix.Data(), run);
102  if (persist) OutputFile = TString::Format("%s_%d_fsim.root",Prefix.Data(), run);
103 
104  gDebug = 0;
105 
106  // choose your event generator
107  Bool_t UseEvtGenDirect = kTRUE;
108  Bool_t UseFtf = kFALSE;
109  Bool_t UseDpm = kFALSE;
110  Bool_t UseBoxGenerator = kFALSE;
111 
112  // use DPM generator; default: inelastic @ pbarmom = mom
113  if (Decfile.BeginsWith("DPM") && !Decfile.EndsWith(".dec"))
114  {
115  UseEvtGenDirect = kFALSE;
116  UseDpm = kTRUE;
117  }
118 
119  // use FTF generator;
120  if (Decfile.BeginsWith("FTF") && !Decfile.EndsWith(".dec"))
121  {
122  UseEvtGenDirect = kFALSE;
123  UseFtf = kTRUE;
124  }
125 
126  // use BOX generator; defaults
127  Double_t BoxMomMin = 0.05; // minimum momentum for box generator
128  Double_t BoxMomMax = Mom; // maximum " "
129  Double_t BoxThtMin = 0. ; // minimum theta for box generator
130  Double_t BoxThtMax = 180.; // maximum " "
131  Double_t BoxPhiMin = 0. ; // minimum phi for box generator
132  Double_t BoxPhiMax = 360.; // maximum " "
133  Bool_t BoxCosTht = false; // isotropic in cos(theta) instead theta
134 
135  Int_t BoxType = 13; // default particle muon
136  Int_t BoxMult = 1; // default particle multiplicity
137  Double_t type=0,mult=0; // ref. parameters for range function
138 
139  if (Decfile.BeginsWith("BOX") && !Decfile.EndsWith(".dec"))
140  {
141  UseEvtGenDirect = kFALSE;
142  UseBoxGenerator = kTRUE;
143  Decfile.ToLower();
144 
145  if (Decfile!="box")
146  {
147  Decfile.ReplaceAll("box","");
148  Decfile.ReplaceAll(" ","");
149  Decfile += ":";
150 
151  while (Decfile.Contains(":"))
152  {
153  TString curpar = Decfile(0,Decfile.Index(":"));
154  Decfile = Decfile(Decfile.Index(":")+1,1000);
155  curpar.ReplaceAll("[","("); curpar.ReplaceAll("]",")");
156 
157  if (curpar.BeginsWith("type(")) {getRange(curpar,type,mult); BoxType = (Int_t)type; BoxMult = (Int_t)mult; }
158  if (curpar.BeginsWith("p(")) getRange(curpar,BoxMomMin,BoxMomMax);
159  if (curpar.BeginsWith("tht(")) getRange(curpar,BoxThtMin,BoxThtMax);
160  if (curpar.BeginsWith("ctht(")) {getRange(curpar,BoxThtMin,BoxThtMax); BoxCosTht=true;}
161  if (curpar.BeginsWith("phi(")) getRange(curpar,BoxPhiMin,BoxPhiMax);
162  }
163  }
164 
165  cout <<"BOX generator range: type["<<BoxType<<","<<BoxMult<<"] p["<<BoxMomMin<<","<<BoxMomMax<<"] tht["<<BoxThtMin<<","<<BoxThtMax<<"]"<<(BoxCosTht?"*":"")<<" phi["<<BoxPhiMin<<","<<BoxPhiMax<<"]"<<endl;
166  }
167 
168 
169  // Start a stop watch
170  TStopwatch timer;
171  timer.Start();
172 
173  // Create the Simulation run manager
174  // --------------------------------
175  FairRunSim *fRun = new FairRunSim();
176  fRun->SetOutputFile(OutputFile.Data());
177  fRun->SetGenerateRunInfo(kFALSE);
178  if (!persist) fRun->SetUserConfig(BaseDir+"/tutorials/analysis/g3ConfigNoMC.C");
179 
180  FairLogger::GetLogger()->SetLogToFile(kFALSE);
181 
182 
183  // Create and Set Event Generator
184  // -------------------------------
186  if (!usePndEventFilter) primGen->SetVerbose(0);
187  //FairPrimaryGenerator* primGen = new FairPrimaryGenerator();
188  fRun->SetGenerator(primGen);
189  fRun->SetName("TGeant3");
190 
191  // Box Generator
192  if(UseBoxGenerator)
193  {
194  PndBoxGenerator* boxGen = new PndBoxGenerator(BoxType, BoxMult);
195  boxGen->SetDebug(0);
196 
197  boxGen->SetPRange(BoxMomMin,BoxMomMax); // GeV/c
198  boxGen->SetPhiRange(BoxPhiMin, BoxPhiMax); // Azimuth angle range [degree]
199  boxGen->SetThetaRange(BoxThtMin, BoxThtMax); // Polar angle in lab system range [degree]
200 
201  if (BoxCosTht) boxGen->SetCosTheta();
202 
203  boxGen->SetXYZ(0., 0., 0.); //cm
204  primGen->AddGenerator(boxGen);
205  }
206 
207  // DPM Generator
208  if(UseDpm)
209  {
210  int mode = 0;
211  if (Decfile=="DPM1") mode = 1;
212  if (Decfile=="DPM2") mode = 2;
213 
214  PndDpmDirect *Dpm= new PndDpmDirect(Mom,mode); // 0 = inelastic, 1 = inelastic & elastic, 2 = elastic
215  // since fastsim doesn't have a transport, let all long-living resonances decay by the generator
216  Dpm->SetUnstable(111); // pi0
217  Dpm->SetUnstable(310); // K_S0
218  Dpm->SetUnstable(311); // K0
219  Dpm->SetUnstable(-311); // K0bar
220  Dpm->SetUnstable(3122); // Lambda0
221  Dpm->SetUnstable(-3122); // anti-Lambda0
222  Dpm->SetUnstable(221); // eta*/
223  Dpm->SetUnstable(3222); // Sigma+
224  Dpm->SetUnstable(-3222); // anti-Sigma-
225  Dpm->SetUnstable(3334); // Omega-
226  primGen->AddGenerator(Dpm);
227  }
228 
229  // FTF Generator
230  if(UseFtf)
231  {
232  bool noelastic = true;
233  if (Decfile=="FTF1") noelastic=false;
234  PndFtfDirect *Ftf = new PndFtfDirect("anti_proton", "G4_H", 1, "ftfp", Mom, 0, noelastic);
235  primGen->AddGenerator(Ftf);
236  }
237 
238  // EvtGen Generator
239  if(UseEvtGenDirect)
240  {
241  TString Resonance=getInitialResonance(Decfile);
242  Resonance.ReplaceAll("pbp","pbarpSystem");
243 
244  PndEvtGenDirect *EvtGen = new PndEvtGenDirect(Resonance, Decfile.Data(), Mom);
245  EvtGen->SetStoreTree(kTRUE);
246  primGen->AddGenerator(EvtGen);
247  }
248 
249  // ------------- switch off the transport of particles
250  primGen->DoTracking(kFALSE);
251 
252  //---------------------Create and Set the Field(s)----------
253  PndMultiField *fField= new PndMultiField("AUTO");
254  fRun->SetField(fField);
255 
256  // Setup the Fast Simulation Task
257  //-----------------------------
258  PndFastSim* fastSim = new PndFastSim(persist);
259 
260  // increasing verbosity increases the amount of console output (mainly for debugging)
261  fastSim->SetVerbosity(0);
262 
263  // set PANDA event filters
264  //-----------------------------
265  if (usePndEventFilter)
266  {
267  cout <<"Using FairEventFilter"<<endl;
268  primGen->SetFilterMaxTries(100000);
269 
271  LambdaFilter->AndMinPdgCodes(1,3122,-3122);
272  //pi0Filter->AndThetaRange(10., 30.);
273  primGen->AndFilter(LambdaFilter);
274 
275  //FairEvtFilterOnSingleParticleCounts* chrgFilter = new FairEvtFilterOnSingleParticleCounts("chrgFilter");
276  //chrgFilter->AndMinCharge(4, FairEvtFilter::kCharged);
277  //primGen->AndFilter(chrgFilter);
278 
279  //FairEvtFilterOnCounts* neutFilter = new FairEvtFilterOnCounts("neutFilter");
280  //neutFilter->AndMaxCharge(4, FairEvtFilter::kNeutral);
281  //primGen->AndFilter(neutFilter);
282  /*
283  FairEvtFilterOnSingleParticleCounts* pi0Filter = new FairEvtFilterOnSingleParticleCounts("pi0Filter");
284  pi0Filter->AndMinPdgCodes(1,111);
285  pi0Filter->AndThetaRange(10., 30.);
286  // require a theta angle of at most 10 degree
287  primGen->AndFilter(pionCand);
288  */
289 
290  /*PndEvtFilterOnInvMassCounts* eeInv= new PndEvtFilterOnInvMassCounts("eeInvMFilter");
291  //eeInv->SetVerbose();//highest commenting level of the FairEvtFilterOnCounts
292  eeInv->SetPdgCodesToCombine( 13, -13);
293  eeInv->SetMinMaxInvMass( 2.5, 3.3 );
294  eeInv->SetMinMaxCounts(1,10000);
295  primGen->AndFilter(eeInv); //add filter to fFilterList
296  */
297 
298  }
299 
300  // set event filters
301  //-----------------------------
302  if (useEventFilter)
303  {
304  // Filters are:
305  // -----------
306  // fastSim->SetMultFilter(type, min, max);
307  // requires min <= mult <= max
308 
309  // available types are:
310 
311  // "+" : positive charged particles
312  // "-" : negative charged particles
313  // "gam" : gammas
314  // "pi0" : pi0 candidates ( -> 2 gammas); mass window 0.135 +- 0.03 GeV
315  // "eta" : eta candidates ( -> 2 gammas); mass window 0.547 +- 0.04 GeV
316  // "ks" : K_S candidates ( -> pi+ pi-); mass window 0.497 +- 0.04 GeV
317 
318  //fastSim->SetMultFilter("ks", 1,1000); // at least 1 KS
319 
320 // fastSim->SetMultFilter("+", 2,1000); // at least 2 trk+
321 // fastSim->SetMultFilter("-", 2,1000); // at least 2 trk-
322 // fastSim->SetMultFilter("gam", 0, 4); // at most 4 gammas
323 
324  // fastSim->SetInvMassFilter(comb, m_min, m_max, mult);
325 
326  // requires at least mult combined candidates with m_min < m < m_max
327 
328  // comb is a TString describing the combinatoric
329  // - particle codes are: e+ e- mu+ mu- pi+ pi- k+ k- p+ p- gam pi0 ks eta
330  // - codes must be separated with a single blank
331  // - for charged final states only the mass is set; no pdg code selection is done!
332  // - optional a 'cc' added at the end of also takes into account charge conjugation
333 
334  // Examples:
335  // - ("k+ k-", 0.98, 1.1, 2) : forms K+ K- candidate and requires >=2 in the given window
336  // - ("ks k+ pi- cc", 2.8, 3.2,1 ) : forms ks k+ pi- / ks k- pi+ cands and req. at least one in window
337 
338  fastSim->SetInvMassFilter("mu+ mu-",2.5,3.3,1); // look for J/psi -> e+ e- candidate
339  }
340 
341  // enable the merging of neutrals if they have similar direction
342  //-----------------------------
343  fastSim->MergeNeutralClusters(mergeNeutrals);
344 
345  // enable bremsstahlung loss for electrons
346  //-----------------------------
347  fastSim->EnableElectronBremsstrahlung(electronBrems);
348 
349  //enable the producting of parametrized neutral (hadronic) split offs
350  // generate electro-magnetic / hadronic split offs in the EMC? switch off when running w/o EMC
351 
352  if (enableSplitoff)
353  fastSim->EnableSplitoffs(splitpars.Data());
354 
355  fastSim->SetUseFlatCov(true);
356  // -----------------------------------------------------------------------------------
357  //Tracking: Set up in parts of theta coverage. All modelled by PndFsmSimpleTracker.
358  // Mind: Numbers on resolution (pRes,thtRes,phiRes) and efficiency are guessed
359  // -----------------------------------------------------------------------------------
360  if (SwMvdGem) // MVD and GEM are enabled; combined tracking available
361  {
362  // - (Full Panda Tracking: STT MVD GEM FTS)
363 
364  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");
365  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");
366  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");
367  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");
368  }
369  else // MVD and GEM are disabled; only STT tracking in central region
370  {
371  // - STT alone:
372  fastSim->AddDetector("ScSttAlone", "thtMin=133.6 thtMax=159.5 ptmin=0.1 pmin=0.0 pRes=0.04 thtRes=0.006 phiRes=0.007 efficiency=0.25");
373  fastSim->AddDetector("ScSttAlone2", "thtMin=20.9 thtMax=133.6 ptmin=0.1 pmin=0.0 pRes=0.04 thtRes=0.006 phiRes=0.007 efficiency=0.80");
374  fastSim->AddDetector("ScSttAlone3", "thtMin=7.8 thtMax=20.9 ptmin=0.1 pmin=0.0 pRes=0.04 thtRes=0.006 phiRes=0.007 efficiency=0.25");
375  }
376 
377  if (SwFwdSpec) // Fwd spectrometer enabled -> use Fwd tracking system
378  {
379  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");
380  }
381 
382  // -----------------------------------------------------------------------------------
383  // Vertexing
384  // -----------------------------------------------------------------------------------
385  if (SwMvdGem) // MVD and GEM are enabled -> better vertexing in central region
386  {
387  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
388  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
389  }
390  else // MVD and GEM are disabled -> no good vertexing at all
391  {
392  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
393  }
394  // -----------------------------------------------------------------------------------
395  // EM Calorimeters w/ default parameters
396  // (don't have to be set, just to list the available parameters
397  // -----------------------------------------------------------------------------------
398 
399  fastSim->AddDetector("EmcFwCap", "thtMin=10.0 thtMax=22.0 Emin=0.01 dist=2.5");
400  fastSim->AddDetector("EmcBwCap", "thtMin=142.0 thtMax=160.0 Emin=0.01 dist=0.7");
401 
402  if (SwEmcBar)
403  {
404  // EmcBarrel also allows to set phiMin and phiMax and can be added multiple times as EmcBarrel1, EmcBarrel2, etc.
405  // Should be made constistent with EmcPidBarrel below
406  fastSim->AddDetector("EmcBarrel","thtMin=22.0 thtMax=142.0 Emin=0.01 barrelRadius=0.5");
407  }
408 
409  if (SwFwdSpec) // Fwd spectrometer enabled -> use Fwd EMC
410  {
411  fastSim->AddDetector("EmcFS", "thtMin=0.05 thtMax=10.0 aPar=0.013 bPar=0.0283 Emin=0.01 dist=8.2");
412  }
413 
414  // -----------------------------------------------------------------------------------
415  // PID
416  // -----------------------------------------------------------------------------------
417 
418  // PID detectors being always in: STT, MUO Barrel, EMC FwdCap, EMC BwdCap
419  //Note: A dEdX parametrization from 2008
420  fastSim->AddDetector("SttPid","thtMin=7.8 thtMax=159.5 ptmin=0.1 dEdxRes=0.15 efficiency=1.");
421  fastSim->AddDetector("ScMdtPidBarrel", "thtMin=10.0 thtMax=130.0 pmin=0.5 efficiency=0.95 misId=0.01");
422  fastSim->AddDetector("ScEmcPidFwCap", "thtMin=10.0 thtMax=22.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
423  fastSim->AddDetector("ScEmcPidBwCap", "thtMin=142.0 thtMax=160.0 ptmin=0.0 pmin=0.0 efficiency=1.0");
424 
425  if (SwMvdGem) // MVD and GEM are enabled -> MVD PID available
426  {
427  //Note: A Bethe-Bloch-Landau-Gauss Prametrization from 2008
428  fastSim->AddDetector("MvdPid","thtMin=5. thtMax=133.6 ptmin=0.1 dEdxResMulti=1. efficiency=1.");
429  }
430 
431  if (SwEmcBar) // EMC Barrel enable -> EMC barrel PID available
432  {
433  fastSim->AddDetector("ScEmcPidBarrel", "thtMin=22.0 thtMax=142.0 ptmin=0.2 pmin=0.0 efficiency=1.0");
434  }
435 
436  if (SwDrc) // Barrel DIRC enabled
437  {
438  fastSim->AddDetector("DrcBarrel","thtMin=22.0 thtMax=140.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
439  }
440 
441  if (SwDsc) // Disc DIRC enabled
442  {
443  fastSim->AddDetector("DrcDisc","thtMin=5.0 thtMax=22.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
444  }
445 
446  if (SwFwdSpec) // Fwd spectrometer enabled -> use RICH, FwdMUO and EMC FS
447  {
448  fastSim->AddDetector("ScEmcPidFS", "thtMin=0.5 thtMax=10.0 ptmin=0.0 pmin=0.5 efficiency=1.0");
449  fastSim->AddDetector("Rich","angleXMax=5.0 angleYMax=10.0 dthtc=0.01 nPhotMin=5 effNPhotons=0.075");
450  fastSim->AddDetector("ScMdtPidForward","thtMin=0.0 thtMax=10.0 pmin=0.5 efficiency=0.95 misId=0.01");
451  }
452 
453 
454  fRun->AddTask(fastSim);
455 
456 
457  // ***********************
458  // *** SoftTriggerTask ***
459  // ***********************
460 
461  if (runST)
462  {
463  // this file contains the trigger line definitions
464  TString triggercfg = TString(gSystem->Getenv("VMCWORKDIR"))+"/softrig/triggerlines_fsim.cfg";
465 
466  PndSoftTriggerTask *stTask = new PndSoftTriggerTask(Mom, 0, 0, triggercfg);
467  stTask->SetFastSimDefaults();
468  stTask->SetQAEvent();
469  fRun->AddTask(stTask);
470  }
471 
472  // ***********************
473  // *** SoftTriggerTask ***
474  // ***********************
475 
476 
477 
478  // *****************************
479  // *** PndSimpleCombinerTask ***
480  // *****************************
481 
482  if (doreco)
483  {
484  PndSimpleCombinerTask *scTask = new PndSimpleCombinerTask(anadecay, anaparms+":algo=PidChargedProbability",Mom, run, runmode);
485  scTask->SetPidAlgo("PidChargedProbability");
486  fRun->AddTask(scTask);
487  }
488 
489  // *****************************
490  // *** PndSimpleCombinerTask ***
491  // *****************************
492 
493 
494 
495  // *****************************
496  // *** PndParticleQATask ***
497  // *****************************
498 
499  if (partQA)
500  {
501  PndParticleQATask *partQaTask = new PndParticleQATask(kTRUE,chrg,neut,mc); // particle QA task for FastSim
502  fRun->AddTask(partQaTask);
503  }
504 
505  // *****************************
506  // *** PndParticleQATask ***
507  // *****************************
508 
509 
510  //------------------------- Initialize the RUN -----------------
511  fRun->Init();
512 
513  //------------------------- Run the Simulation -----------------
514  fRun->Run(nEvents);
515 
516  //------------------------- Write Filter Info to File -----------
517  if (usePndEventFilter) primGen->WriteEvtFilterStatsToRootFile();
518 
519  //------------------------Print some info and exit----------------
520  timer.Stop();
521  Double_t rtime = timer.RealTime();
522  Double_t ctime = timer.CpuTime();
523  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
524 }
525 
526 void getRange(TString par, double &min, double &max)
527 {
528  par.ReplaceAll(" ","");
529  par = par(par.Index("(")+1, par.Length()-par.Index("(")-2);
530 
531  TString smin=par, smax=par;
532 
533  if (par.Contains(","))
534  {
535  smin = par(0,par.Index(","));
536  smax = par(par.Index(",")+1,1000);
537  }
538 
539  min = smin.Atof();
540  max = smax.Atof();
541 
542  //if (min>max) {double tmp=min;min=max;max=tmp;}
543 }
544 
545 TString getInitialResonance(TString &fEvtGenFile)
546 {
547 
548  TString IniRes="";
549 
550  if (fEvtGenFile.Contains(":")) // is the initial resonance provide as <decfile>.dec:iniRes ?
551  {
552  IniRes = fEvtGenFile(fEvtGenFile.Index(":")+1,1000);
553  fEvtGenFile = fEvtGenFile(0,fEvtGenFile.Index(":"));
554  }
555 
556  if (IniRes=="") // we need to search the decay file
557  {
558  std::ifstream fs(fEvtGenFile.Data());
559  char line[250];
560 
561  while (fs)
562  {
563  fs.getline(line,249);
564  TString s(line);
565  s.ReplaceAll("\r","");
566  if (IniRes=="" && s.Contains("Decay "))
567  {
568  if (s.Contains("#")) s=s(0,s.Index("#"));
569  s.ReplaceAll("Decay ","");
570  s.ReplaceAll(" ","");
571  IniRes = s;
572  }
573  }
574  fs.close();
575  }
576 
577  return IniRes;
578 }
PndMultiField * fField
Definition: sim_emc_apd.C:97
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
void SetThetaRange(Double32_t thetamin=0, Double32_t thetamax=90)
Int_t run
Definition: autocutx.C:47
void WriteEvtFilterStatsToRootFile(TFile *outputFile=NULL)
Writes all relevant event filter information to the output root file.
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void AndFilter(FairEvtFilter *filter)
Register a non-veto event filter using a logical AND to connect with previously defined non-veto even...
bool AddDetector(std::string name, std::string params="")
Definition: PndFastSim.cxx:313
TLorentzVector s
Definition: Pnd2DStar.C:50
void SetXYZ(Double32_t x=0, Double32_t y=0, Double32_t z=0)
void SetPidAlgo(TString algo)
Double_t par[3]
TString getInitialResonance(TString &fEvtGenFile)
void SetFilterMaxTries(Int_t maxTries=99999)
Define the maximum number of times that this object should try to find an event which suits all event...
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
Bool_t AndMinPdgCodes(Int_t min, Int_t pdgCode1, Int_t pdgCode2=kInvalidPdgCode, Int_t pdgCode3=kInvalidPdgCode, Int_t pdgCode4=kInvalidPdgCode, Int_t pdgCode5=kInvalidPdgCode, Int_t pdgCode6=kInvalidPdgCode, Int_t pdgCode7=kInvalidPdgCode, Int_t pdgCode8=kInvalidPdgCode)
void getRange(TString par, double &min, double &max)
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)
Primary generator with added event filtering capabilities.
bool MergeNeutralClusters(bool merge=true, double par=0.389)
Definition: PndFastSim.h:61
void SetInvMassFilter(TString filter, double min, double max, int mult=1)
Definition: PndFastSim.cxx:373
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
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
void SetDebug(Bool_t debug=0)
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
double X
Definition: anaLmdDigi.C:68
void SetVerbosity(int vb)
Definition: PndFastSim.h:59
void quickfsimana(TString Prefix="", TString Decfile="", Float_t Mom=0., TString anadecay="", Int_t nEvents=1000, TString anaparms="", bool runST=false, int run=0, int runmode=0)
void SetPRange(Double32_t pmin=0, Double32_t pmax=10)
TLorentzVector fIni
Definition: full_core_ntp.C:29
Double_t rtime
Definition: hit_dirc.C:113
Double_t mult
void SetPhiRange(Double32_t phimin=0, Double32_t phimax=360)
void SetVerbose(Int_t verbose=12)
Set the level of commenting output.