FairRoot/PandaRoot
outdated/run/ana_complete.C
Go to the documentation of this file.
1 class RhoCandList;
2 class RhoCandidate;
5 class PndAnalysis;
6 
7 // *** routine to only keep PID matched candidates in list
9 {
10  int removed = 0;
11 
12  for (int ii=l.GetLength()-1;ii>=0;--ii)
13  {
14  if ( !(ana->McTruthMatch(l[ii])) )
15  {
16  l.Remove(l[ii]);
17  removed++;
18  }
19  }
20 
21  return removed;
22 }
23 
25 {
26  TLorentzVector lv=c->P4();
27 
28  cout <<c->PdgCode()<<" ("<<lv.X()<<"/"<<lv.Y()<<"/"<<lv.Z()<<"/"<<lv.E()<<")"<<endl;
29 }
30 
31 void countDoubles(RhoCandList &l, int &n1, int &n2, int &n3)
32 {
33  int n_smc = 0;
34  int n_strk = 0;
35  int n_both = 0;
36  double d = 0.00001;
37 
38  for (int i=0;i<l.GetLength()-1;++i)
39  {
40  for (int j=i+1;j<l.GetLength();++j)
41  {
42  TLorentzVector dl = l[i]->P4() - l[j]->P4();
43 
44  bool chkmc = (l[i]->GetMcTruth()==l[j]->GetMcTruth());
45  bool chktrk = (fabs(dl.X())<d) && (fabs(dl.Y())<d) && (fabs(dl.Z())<d) && (fabs(dl.E())<d);
46  if (chkmc) n_smc++;
47  if (chktrk) n_strk++;
48  if (chktrk && chkmc) n_both++;
49  }
50  }
51  n1 = n_strk;
52  n2 = n_smc;
53  n3 = n_both;
54 }
55 
56 int ana_complete(int nevts=0)
57 {
58  TDatabasePDG::Instance()->AddParticle("pbarpSystem","pbarpSystem",1.9,kFALSE,0.1,0,"",88888);
59  TStopwatch timer;
60  // *** some variables
61  int i=0,j=0, k=0, l=0;
62  gStyle->SetOptFit(1011);
63 
64  // *** the output file for FairRunAna
65  TString OutFile="output.root";
66 
67  // *** the files coming from the simulation
68  TString inPidFile = "pid_complete.root"; // this file contains the PndPidCandidates and McTruth
69  TString inParFile = "simparams.root";
70 
71  // *** PID table with selection thresholds; can be modified by the user
72  TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
73 
74  // *** initialization
75  FairLogger::GetLogger()->SetLogToFile(kFALSE);
76  FairRunAna* fRun = new FairRunAna();
77  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
78  fRun->SetInputFile(inPidFile);
79 
80  // *** setup parameter database
81  FairParRootFileIo* parIO = new FairParRootFileIo();
82  parIO->open(inParFile);
83  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
84  parIOPid->open(pidParFile.Data(),"in");
85 
86  rtdb->setFirstInput(parIO);
87  rtdb->setSecondInput(parIOPid);
88  rtdb->setOutput(parIO);
89 
90  fRun->SetOutputFile(OutFile);
91  fRun->Init();
92 
93  // *** create an output file for all histograms
94  TFile *out = TFile::Open("output_ana.root","RECREATE");
95 
96  // *** create some histograms
97  TH1F *hmomtrk = new TH1F("hmomtrk","track momentum (all)",200,0,5);
98  TH1F *hthttrk = new TH1F("hthttrk","track theta (all)",200,0,3.1415);
99 
100  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
101  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
102 
103  TH1F *hjpsim_lpid = new TH1F("hjpsim_lpid","J/#psi mass (loose pid)",200,0,4.5);
104  TH1F *hpsim_lpid = new TH1F("hpsim_lpid","#psi(2S) mass (loose pid)",200,0,5);
105 
106  TH1F *hjpsim_tpid = new TH1F("hjpsim_tpid","J/#psi mass (tight pid)",200,0,4.5);
107  TH1F *hpsim_tpid = new TH1F("hpsim_tpid","#psi(2S) mass (tight pid)",200,0,5);
108 
109  TH1F *hjpsim_trpid = new TH1F("hjpsim_trpid","J/#psi mass (true pid)",200,0,4.5);
110  TH1F *hpsim_trpid = new TH1F("hpsim_trpid","#psi(2S) mass (true pid)",200,0,5);
111 
112 
113  TH1F *hjpsim_ftm = new TH1F("hjpsim_ftm","J/#psi mass (full truth match)",200,0,4.5);
114  TH1F *hpsim_ftm = new TH1F("hpsim_ftm","#psi(2S) mass (full truth match)",200,0,5);
115 
116  TH1F *hjpsim_nm = new TH1F("hjpsim_nm","J/#psi mass (no truth match)",200,0,4.5);
117  TH1F *hpsim_nm = new TH1F("hpsim_nm","#psi(2S) mass (no truth match)",200,0,5);
118 
119  TH1F *hjpsim_diff = new TH1F("hjpsim_diff","J/#psi mass diff to truth",100,-2,2);
120  TH1F *hpsim_diff = new TH1F("hpsim_diff","#psi(2S) mass diff to truth",100,-2,2);
121 
122 
123  TH1F *hjpsim_vf = new TH1F("hjpsim_vf","J/#psi mass (vertex fit)",200,0,4.5);
124  TH1F *hjpsim_4cf = new TH1F("hjpsim_4cf","J/#psi mass (4C fit)",200,0,4.5);
125  TH1F *hjpsim_mcf = new TH1F("hjpsim_mcf","J/#psi mass (mass constraint fit)",200,0,4.5);
126 
127  TH1F *hjpsi_chi2_vf = new TH1F("hjpsi_chi2_vf", "J/#psi: #chi^{2} vertex fit",100,0,10);
128  TH1F *hpsi_chi2_4c = new TH1F("hpsi_chi2_4c", "#psi(2S): #chi^{2} 4C fit",100,0,250);
129  TH1F *hjpsi_chi2_mf = new TH1F("hjpsi_chi2_mf", "J/#psi: #chi^{2} mass fit",100,0,10);
130 
131  TH1F *hjpsi_prob_vf = new TH1F("hjpsi_prob_vf", "J/#psi: Prob vertex fit",100,0,1);
132  TH1F *hpsi_prob_4c = new TH1F("hpsi_prob_4c", "#psi(2S): Prob 4C fit",100,0,1);
133  TH1F *hjpsi_prob_mf = new TH1F("hjpsi_prob_mf", "J/#psi: Prob mass fit",100,0,1);
134 
135  TH2F *hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-2,2,100,-2,2);
136 
137  //
138  // Now the analysis stuff comes...
139  //
140 
141 
142  // *** the data reader object
143  PndAnalysis* theAnalysis = new PndAnalysis();
144  if (nevts==0) nevts= theAnalysis->GetEntries();
145 
146  // *** RhoCandLists for the analysis
147  RhoCandList chrg, muplus, muminus, piplus, piminus, jpsi, psi2s;
148 
149  // *** Mass selector for the jpsi cands
150  double m0_jpsi = TDatabasePDG::Instance()->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
151  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",m0_jpsi,1.0);
152 
153  // *** the lorentz vector of the initial psi(2S)
154  TLorentzVector ini(0, 0, 6.231552, 7.240065);
155 
156  // ***
157  // the event loop
158  // ***
159 
160  int cntdbltrk=0, cntdblmc=0, cntdblboth=0, cnttrk=0, cnt_dbl_jpsi=0, cnt_dbl_psip=0;
161 
162  while (theAnalysis->GetEvent() && i++<nevts)
163  {
164  if ((i%100)==0) cout<<"evt " << i << endl;
165 
166  // *** Select with no PID info ('All'); type and mass are set
167  theAnalysis->FillList(chrg, "Charged");
168  theAnalysis->FillList(muplus, "MuonAllPlus");
169  theAnalysis->FillList(muminus, "MuonAllMinus");
170  theAnalysis->FillList(piplus, "PionAllPlus");
171  theAnalysis->FillList(piminus, "PionAllMinus");
172 
173  // *** momentum and theta histograms
174  for (j=0;j<muplus.GetLength();++j)
175  {
176  hmomtrk->Fill(muplus[j]->P());
177  hthttrk->Fill(muplus[j]->P4().Theta());
178  }
179  for (j=0;j<muminus.GetLength();++j)
180  {
181  hmomtrk->Fill(muminus[j]->P());
182  hthttrk->Fill(muminus[j]->P4().Theta());
183  }
184 
185  cnttrk += chrg.GetLength();
186 
187  int n1, n2, n3;
188 
189  countDoubles(chrg,n1,n2,n3);
190  cntdbltrk += n1;
191  cntdblmc += n2;
192  cntdblboth += n3;
193 
194  // *** combinatorics for J/psi -> mu+ mu-
195  jpsi.Combine(muplus, muminus);
196 
197 
198  // ***
199  // *** do the TRUTH MATCH for jpsi
200  // ***
201  jpsi.SetType(443);
202 
203  int nm = 0;
204  for (j=0;j<jpsi.GetLength();++j)
205  {
206  hjpsim_all->Fill( jpsi[j]->M() );
207 
208  if (theAnalysis->McTruthMatch(jpsi[j]))
209  {
210  nm++;
211  hjpsim_ftm->Fill( jpsi[j]->M() );
212  hjpsim_diff->Fill( jpsi[j]->GetMcTruth()->M() - jpsi[j]->M() );
213  }
214  else
215  hjpsim_nm->Fill( jpsi[j]->M() );
216  }
217 
218  if (nm>1) cnt_dbl_jpsi++;
219  // ***
220  // *** do VERTEX FIT (J/psi)
221  // ***
222  for (j=0;j<jpsi.GetLength();++j)
223  {
224  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
225  vtxfitter.Fit();
226 
227  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
228  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
229  hjpsi_chi2_vf->Fill(chi2_vtx);
230  hjpsi_prob_vf->Fill(prob_vtx);
231 
232  if ( prob_vtx > 0.01 ) // when good enough, fill some histos
233  {
234  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
235  TVector3 jVtx=jfit->Pos(); // and the decay vertex position
236 
237  hjpsim_vf->Fill(jfit->M());
238  hvpos->Fill(jVtx.X(),jVtx.Y());
239  }
240  }
241 
242  // *** some rough mass selection
243  jpsi.Select(jpsiMassSel);
244 
245  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
246  psi2s.Combine(jpsi, piplus, piminus);
247 
248 
249  // ***
250  // *** do the TRUTH MATCH for psi(2S)
251  // ***
252  psi2s.SetType(88888);
253 
254  nm = 0;
255  for (j=0;j<psi2s.GetLength();++j)
256  {
257  hpsim_all->Fill( psi2s[j]->M() );
258 
259  if (theAnalysis->McTruthMatch(psi2s[j]))
260  {
261  nm++;
262  hpsim_ftm->Fill( psi2s[j]->M() );
263  hpsim_diff->Fill( psi2s[j]->GetMcTruth()->M() - psi2s[j]->M() );
264  }
265  else
266  hpsim_nm->Fill( psi2s[j]->M() );
267  }
268  if (nm>1) cnt_dbl_psip++;
269 
270 
271  // ***
272  // *** do 4C FIT (initial psi(2S) system)
273  // ***
274  for (j=0;j<psi2s.GetLength();++j)
275  {
276  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
277  fitter.Add4MomConstraint(ini); // set 4 constraint
278  fitter.Fit(); // do fit
279 
280  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
281  double prob_4c = fitter.GetProb(); // access probability of fit
282  hpsi_chi2_4c->Fill(chi2_4c);
283  hpsi_prob_4c->Fill(prob_4c);
284 
285  if ( prob_4c > 0.01 ) // when good enough, fill some histo
286  {
287  RhoCandidate *jfit = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
288 
289  hjpsim_4cf->Fill(jfit->M());
290  }
291  }
292 
293 
294  // ***
295  // *** do MASS CONSTRAINT FIT (J/psi)
296  // ***
297  for (j=0;j<jpsi.GetLength();++j)
298  {
299  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
300  mfitter.AddMassConstraint(m0_jpsi); // add the mass constraint
301  mfitter.Fit(); // do fit
302 
303  double chi2_m = mfitter.GetChi2(); // get chi2 of fit
304  double prob_m = mfitter.GetProb(); // access probability of fit
305  hjpsi_chi2_mf->Fill(chi2_m);
306  hjpsi_prob_mf->Fill(prob_m);
307 
308  if ( prob_m > 0.01 ) // when good enough, fill some histo
309  {
310  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
311  hjpsim_mcf->Fill(jfit->M());
312  }
313  }
314 
315 
316  // ***
317  // *** TRUE PID combinatorics
318  // ***
319 
320  // *** do MC truth match for PID type
321  SelectTruePid(theAnalysis, muplus);
322  SelectTruePid(theAnalysis, muminus);
323  SelectTruePid(theAnalysis, piplus);
324  SelectTruePid(theAnalysis, piminus);
325 
326  // *** all combinatorics again with true PID
327  jpsi.Combine(muplus, muminus);
328  for (j=0;j<jpsi.GetLength();++j) hjpsim_trpid->Fill( jpsi[j]->M() );
329  jpsi.Select(jpsiMassSel);
330 
331  psi2s.Combine(jpsi, piplus, piminus);
332  for (j=0;j<psi2s.GetLength();++j) hpsim_trpid->Fill( psi2s[j]->M() );
333 
334 
335  // ***
336  // *** LOOSE PID combinatorics
337  // ***
338 
339  // *** and again with PidAlgoMvd;PidAlgoStt;PidAlgoDrc and loose selection
340  theAnalysis->FillList(muplus, "MuonLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc;PidAlgoDisc;PidAlgoMdtHardCuts");
341  theAnalysis->FillList(muminus, "MuonLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc;PidAlgoDisc;PidAlgoMdtHardCuts");
342  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc;PidAlgoDisc");
343  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc;PidAlgoDisc");
344 
345  jpsi.Combine(muplus, muminus);
346  for (j=0;j<jpsi.GetLength();++j) hjpsim_lpid->Fill( jpsi[j]->M() );
347  jpsi.Select(jpsiMassSel);
348 
349  psi2s.Combine(jpsi, piplus, piminus);
350  for (j=0;j<psi2s.GetLength();++j) hpsim_lpid->Fill( psi2s[j]->M() );
351 
352 
353  // ***
354  // *** TIGHT PID combinatorics
355  // ***
356 
357  // *** and again with PidAlgoMvd;PidAlgoStt and tight selection
358  theAnalysis->FillList(muplus, "MuonTightPlus", "PidAlgoMdtHardCuts");
359  theAnalysis->FillList(muminus, "MuonTightMinus", "PidAlgoMdtHardCuts");
360  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
361  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
362 
363  jpsi.Combine(muplus, muminus);
364  for (j=0;j<jpsi.GetLength();++j) hjpsim_tpid->Fill( jpsi[j]->M() );
365  jpsi.Select(jpsiMassSel);
366 
367  psi2s.Combine(jpsi, piplus, piminus);
368  for (j=0;j<psi2s.GetLength();++j) hpsim_tpid->Fill( psi2s[j]->M() );
369 
370  }
371 
372  // *** write out all the histos
373  out->cd();
374 
375  hmomtrk->Write();
376  hthttrk->Write();
377 
378  hjpsim_all->Write();
379  hpsim_all->Write();
380  hjpsim_lpid->Write();
381  hpsim_lpid->Write();
382  hjpsim_tpid->Write();
383  hpsim_tpid->Write();
384  hjpsim_trpid->Write();
385  hpsim_trpid->Write();
386 
387  hjpsim_ftm->Write();
388  hpsim_ftm->Write();
389  hjpsim_nm->Write();
390  hpsim_nm->Write();
391 
392  hpsim_diff->Write();
393  hjpsim_diff->Write();
394 
395  hjpsim_vf->Write();
396  hjpsim_4cf->Write();
397  hjpsim_mcf->Write();
398 
399  hjpsi_chi2_vf->Write();
400  hpsi_chi2_4c->Write();
401  hjpsi_chi2_mf->Write();
402 
403  hjpsi_prob_vf->Write();
404  hpsi_prob_4c->Write();
405  hjpsi_prob_mf->Write();
406 
407  hvpos->Write();
408 
409  out->Save();
410 
411  timer.Stop();
412  Double_t rtime = timer.RealTime();
413  Double_t ctime = timer.CpuTime();
414  cout << endl << endl;
415  cout << "Macro finished successfully." << endl;
416  cout << "Real time " << rtime << " s, CPU time " << ctime << " s" << endl;
417  cout << endl;
418  // ------------------------------------------------------------------------
419  cout << " Test passed" << endl;
420  cout << " All ok " << endl;
421  return 0;
422 
423 }
Int_t GetEntries()
int SelectTruePid(PndAnalysis *ana, RhoCandList &l)
void AddMassConstraint(double mass)
TObjArray * d
Int_t i
Definition: run_full.C:25
TVector3 Pos() const
Definition: RhoCandidate.h:186
Int_t GetLength() const
Definition: RhoCandList.h:46
TH2F * hvpos
Definition: plot_eta_c.C:47
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
CandList piplus
void Add4MomConstraint(TLorentzVector lv)
FairRunAna * fRun
Definition: hit_dirc.C:58
CandList muplus
void Combine(RhoCandList &l1, RhoCandList &l2)
Double_t
PndAnalysis(TString tname1="", TString tname2="", TString algnamec="PidAlgoIdealCharged", TString algnamen="PidAlgoIdealNeutral")
Definition: PndAnalysis.cxx:48
TStopwatch timer
Definition: hit_dirc.C:51
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
TLorentzVector P4() const
Definition: RhoCandidate.h:195
double GetProb() const
Definition: RhoFitterBase.h:50
void SetType(const TParticlePDG *pdt, int start=0)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TFile * out
Definition: reco_muo.C:20
Double_t M() const
Int_t Remove(RhoCandidate *)
void countDoubles(RhoCandList &l, int &n1, int &n2, int &n3)
Double_t ctime
Definition: hit_dirc.C:114
GeV c P
CandList piminus
void printCand(TLorentzVector l, TVector3 p)
double GetChi2() const
Definition: RhoFitterBase.h:48
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
Int_t GetEvent(Int_t n=-1)
Bool_t Fit()
CandList muminus
Double_t rtime
Definition: hit_dirc.C:113
int ana_complete(int nevts=0, TString prefix="evtcomplete")