FairRoot/PandaRoot
rho/tut_ana_fast.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 
24 
25 void tut_ana_fast(int nevts = 0, TString prefix = "signal")
26 {
27  // *** some variables
28  int i=0,j=0, k=0, l=0;
29  gStyle->SetOptFit(1011);
30 
31  // *** the output file for FairRunAna
32  TString OutFile="out_dummy.root";
33 
34  // *** the files coming from the simulation
35  TString inPidFile = prefix+"_fast.root"; // this file contains the PndPidCandidates and McTruth
36 
37  // *** PID table with selection thresholds; can be modified by the user
38  TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
39 
40  // *** initialization
41  //FairLogger::GetLogger()->SetLogToFile(kFALSE);
42  FairRunAna* fRun = new FairRunAna();
43  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
44  fRun->SetSource(new FairFileSource(inPidFile));
45  fRun->SetOutputFile(OutFile);
46 
47  // *** setup parameter database
48  FairParRootFileIo* parIO = new FairParRootFileIo();
49  //parIO->open(inParFile);
50  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
51  parIOPid->open(pidParFile.Data(),"in");
52 
53  rtdb->setFirstInput(parIO);
54  rtdb->setSecondInput(parIOPid);
55  rtdb->setOutput(parIO);
56 
57  fRun->SetOutputFile(OutFile);
58  // *** take constant field; needed for PocaVtx
60 
61  fRun->Init();
62 
63  // *** create an output file for all histograms
64  TFile *out = TFile::Open(prefix+"_ana_fast.root","RECREATE");
65 
66  // *** create some histograms
67  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
68  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
69 
70  TH1F *hjpsim_lpid = new TH1F("hjpsim_lpid","J/#psi mass (loose pid)",200,0,4.5);
71  TH1F *hpsim_lpid = new TH1F("hpsim_lpid","#psi(2S) mass (loose pid)",200,0,5);
72 
73  TH1F *hjpsim_tpid = new TH1F("hjpsim_tpid","J/#psi mass (tight pid)",200,0,4.5);
74  TH1F *hpsim_tpid = new TH1F("hpsim_tpid","#psi(2S) mass (tight pid)",200,0,5);
75 
76  TH1F *hjpsim_trpid = new TH1F("hjpsim_trpid","J/#psi mass (true pid)",200,0,4.5);
77  TH1F *hpsim_trpid = new TH1F("hpsim_trpid","#psi(2S) mass (true pid)",200,0,5);
78 
79 
80  TH1F *hjpsim_ftm = new TH1F("hjpsim_ftm","J/#psi mass (full truth match)",200,0,4.5);
81  TH1F *hpsim_ftm = new TH1F("hpsim_ftm","#psi(2S) mass (full truth match)",200,0,5);
82 
83  TH1F *hjpsim_nm = new TH1F("hjpsim_nm","J/#psi mass (no truth match)",200,0,4.5);
84  TH1F *hpsim_nm = new TH1F("hpsim_nm","#psi(2S) mass (no truth match)",200,0,5);
85 
86  TH1F *hjpsim_diff = new TH1F("hjpsim_diff","J/#psi mass diff to truth",100,-2,2);
87  TH1F *hpsim_diff = new TH1F("hpsim_diff","#psi(2S) mass diff to truth",100,-2,2);
88 
89 
90  TH1F *hjpsim_vf = new TH1F("hjpsim_vf","J/#psi mass (vertex fit)",200,0,4.5);
91  TH1F *hjpsim_4cf = new TH1F("hjpsim_4cf","J/#psi mass (4C fit)",200,0,4.5);
92  TH1F *hjpsim_mcf = new TH1F("hjpsim_mcf","J/#psi mass (mass constraint fit)",200,0,4.5);
93 
94  TH1F *hjpsi_chi2_vf = new TH1F("hjpsi_chi2_vf", "J/#psi: #chi^{2} vertex fit",100,0,10);
95  TH1F *hpsi_chi2_4c = new TH1F("hpsi_chi2_4c", "#psi(2S): #chi^{2} 4C fit",100,0,250);
96  TH1F *hjpsi_chi2_mf = new TH1F("hjpsi_chi2_mf", "J/#psi: #chi^{2} mass fit",100,0,10);
97 
98  TH1F *hjpsi_prob_vf = new TH1F("hjpsi_prob_vf", "J/#psi: Prob vertex fit",100,0,1);
99  TH1F *hpsi_prob_4c = new TH1F("hpsi_prob_4c", "#psi(2S): Prob 4C fit",100,0,1);
100  TH1F *hjpsi_prob_mf = new TH1F("hjpsi_prob_mf", "J/#psi: Prob mass fit",100,0,1);
101 
102  TH2F *hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-2,2,100,-2,2);
103 
104  //
105  // Now the analysis stuff comes...
106  //
107 
108 
109  // *** the data reader object
110  PndAnalysis* theAnalysis = new PndAnalysis();
111  if (nevts==0) nevts= theAnalysis->GetEntries();
112 
113  // *** RhoCandLists for the analysis
114  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
115 
116  // *** Mass selector for the jpsi cands
117  double m0_jpsi = TDatabasePDG::Instance()->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
118  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",m0_jpsi,1.0);
119 
120  // *** the lorentz vector of the initial psi(2S)
121  TLorentzVector ini(0, 0, 6.231552, 7.240065);
122 
123  TString pidalgo = "PidChargedProbability";
124 
125  // ***
126  // the event loop
127  // ***
128  while (theAnalysis->GetEvent() && i++<nevts)
129  {
130  if ((i%100)==0) cout<<"evt " << i << endl;
131 
132  // *** Select with no PID info ('All'); type and mass are set
133  theAnalysis->FillList(muplus, "MuonAllPlus" , pidalgo );
134  theAnalysis->FillList(muminus, "MuonAllMinus" , pidalgo );
135  theAnalysis->FillList(piplus, "PionAllPlus" , pidalgo );
136  theAnalysis->FillList(piminus, "PionAllMinus" , pidalgo );
137 
138  // *** combinatorics for J/psi -> mu+ mu-
139  jpsi.Combine(muplus, muminus);
140 
141 
142  // ***
143  // *** do the TRUTH MATCH for jpsi
144  // ***
145  jpsi.SetType(443);
146 
147  for (j=0;j<jpsi.GetLength();++j)
148  {
149  hjpsim_all->Fill( jpsi[j]->M() );
150 
151  if (theAnalysis->McTruthMatch(jpsi[j]))
152  {
153  hjpsim_ftm->Fill( jpsi[j]->M() );
154  hjpsim_diff->Fill( jpsi[j]->GetMcTruth()->M() - jpsi[j]->M() );
155  }
156  else
157  hjpsim_nm->Fill( jpsi[j]->M() );
158  }
159 
160  // ***
161  // *** do VERTEX FIT (J/psi)
162  // ***
163  for (j=0;j<jpsi.GetLength();++j)
164  {
165  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
166  vtxfitter.Fit();
167 
168  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
169  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
170  hjpsi_chi2_vf->Fill(chi2_vtx);
171  hjpsi_prob_vf->Fill(prob_vtx);
172 
173  if ( prob_vtx > 0.01 ) // when good enough, fill some histos
174  {
175  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
176  TVector3 jVtx=jfit->Pos(); // and the decay vertex position
177 
178  hjpsim_vf->Fill(jfit->M());
179  hvpos->Fill(jVtx.X(),jVtx.Y());
180  }
181  }
182 
183  // *** some rough mass selection
184  jpsi.Select(jpsiMassSel);
185 
186  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
187  psi2s.Combine(jpsi, piplus, piminus);
188 
189 
190  // ***
191  // *** do the TRUTH MATCH for psi(2S)
192  // ***
193  psi2s.SetType(88888);
194 
195  for (j=0;j<psi2s.GetLength();++j)
196  {
197  hpsim_all->Fill( psi2s[j]->M() );
198 
199  if (theAnalysis->McTruthMatch(psi2s[j]))
200  {
201  hpsim_ftm->Fill( psi2s[j]->M() );
202  hpsim_diff->Fill( psi2s[j]->GetMcTruth()->M() - psi2s[j]->M() );
203  }
204  else
205  hpsim_nm->Fill( psi2s[j]->M() );
206  }
207 
208 
209  // ***
210  // *** do 4C FIT (initial psi(2S) system)
211  // ***
212  for (j=0;j<psi2s.GetLength();++j)
213  {
214  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
215  fitter.Add4MomConstraint(ini); // set 4 constraint
216  fitter.Fit(); // do fit
217 
218  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
219  double prob_4c = fitter.GetProb(); // access probability of fit
220  hpsi_chi2_4c->Fill(chi2_4c);
221  hpsi_prob_4c->Fill(prob_4c);
222 
223  if ( prob_4c > 0.01 ) // when good enough, fill some histo
224  {
225  RhoCandidate *jfit = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
226 
227  hjpsim_4cf->Fill(jfit->M());
228  }
229  }
230 
231 
232  // ***
233  // *** do MASS CONSTRAINT FIT (J/psi)
234  // ***
235  for (j=0;j<jpsi.GetLength();++j)
236  {
237  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
238  mfitter.AddMassConstraint(m0_jpsi); // add the mass constraint
239  mfitter.Fit(); // do fit
240 
241  double chi2_m = mfitter.GetChi2(); // get chi2 of fit
242  double prob_m = mfitter.GetProb(); // access probability of fit
243  hjpsi_chi2_mf->Fill(chi2_m);
244  hjpsi_prob_mf->Fill(prob_m);
245 
246  if ( prob_m > 0.01 ) // when good enough, fill some histo
247  {
248  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
249  hjpsim_mcf->Fill(jfit->M());
250  }
251  }
252 
253 
254  // ***
255  // *** TRUE PID combinatorics
256  // ***
257 
258  // *** do MC truth match for PID type
259  SelectTruePid(theAnalysis, muplus);
260  SelectTruePid(theAnalysis, muminus);
261  SelectTruePid(theAnalysis, piplus);
262  SelectTruePid(theAnalysis, piminus);
263 
264  // *** all combinatorics again with true PID
265  jpsi.Combine(muplus, muminus);
266  for (j=0;j<jpsi.GetLength();++j) hjpsim_trpid->Fill( jpsi[j]->M() );
267  jpsi.Select(jpsiMassSel);
268 
269  psi2s.Combine(jpsi, piplus, piminus);
270  for (j=0;j<psi2s.GetLength();++j) hpsim_trpid->Fill( psi2s[j]->M() );
271 
272 
273  // ***
274  // *** LOOSE PID combinatorics
275  // ***
276  pidalgo = "MvdPidProbability;SttPidProbability;DrcBarrelProbability";
277  // *** and again with Mvd, Stt, Drc and loose selection
278  theAnalysis->FillList(muplus, "MuonLoosePlus", pidalgo );
279  theAnalysis->FillList(muminus, "MuonLooseMinus", pidalgo );
280  theAnalysis->FillList(piplus, "PionLoosePlus", pidalgo );
281  theAnalysis->FillList(piminus, "PionLooseMinus", pidalgo);
282 
283  jpsi.Combine(muplus, muminus);
284  for (j=0;j<jpsi.GetLength();++j) hjpsim_lpid->Fill( jpsi[j]->M() );
285  jpsi.Select(jpsiMassSel);
286 
287  psi2s.Combine(jpsi, piplus, piminus);
288  for (j=0;j<psi2s.GetLength();++j) hpsim_lpid->Fill( psi2s[j]->M() );
289 
290 
291  // ***
292  // *** TIGHT PID combinatorics
293  // ***
294  TString pidalgomdt = "ScMdtPidBarrelProbability;ScMdtPidForwardProbability";
295  // *** and again with Mvd, Stt and tight selection
296  theAnalysis->FillList(muplus, "MuonTightPlus", pidalgomdt);
297  theAnalysis->FillList(muminus, "MuonTightMinus", pidalgomdt);
298  theAnalysis->FillList(piplus, "PionLoosePlus", pidalgo);
299  theAnalysis->FillList(piminus, "PionLooseMinus", pidalgo);
300 
301  jpsi.Combine(muplus, muminus);
302  for (j=0;j<jpsi.GetLength();++j) hjpsim_tpid->Fill( jpsi[j]->M() );
303  jpsi.Select(jpsiMassSel);
304 
305  psi2s.Combine(jpsi, piplus, piminus);
306  for (j=0;j<psi2s.GetLength();++j) hpsim_tpid->Fill( psi2s[j]->M() );
307 
308  }
309 
310  // *** write out all the histos
311  out->cd();
312 
313  hjpsim_all->Write();
314  hpsim_all->Write();
315  hjpsim_lpid->Write();
316  hpsim_lpid->Write();
317  hjpsim_tpid->Write();
318  hpsim_tpid->Write();
319  hjpsim_trpid->Write();
320  hpsim_trpid->Write();
321 
322  hjpsim_ftm->Write();
323  hpsim_ftm->Write();
324  hjpsim_nm->Write();
325  hpsim_nm->Write();
326 
327  hpsim_diff->Write();
328  hjpsim_diff->Write();
329 
330  hjpsim_vf->Write();
331  hjpsim_4cf->Write();
332  hjpsim_mcf->Write();
333 
334  hjpsi_chi2_vf->Write();
335  hpsi_chi2_4c->Write();
336  hjpsi_chi2_mf->Write();
337 
338  hjpsi_prob_vf->Write();
339  hpsi_prob_4c->Write();
340  hjpsi_prob_mf->Write();
341 
342  hvpos->Write();
343 
344  out->Save();
345 
346 }
Int_t GetEntries()
int SelectTruePid(PndAnalysis *ana, RhoCandList &l)
void AddMassConstraint(double mass)
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)
static void ForceConstantBz(Double_t bz=0.)
Force a constant B field value for all positions.
FairRunAna * fRun
Definition: hit_dirc.C:58
CandList muplus
void Combine(RhoCandList &l1, RhoCandList &l2)
PndAnalysis(TString tname1="", TString tname2="", TString algnamec="PidAlgoIdealCharged", TString algnamen="PidAlgoIdealNeutral")
Definition: PndAnalysis.cxx:48
void tut_ana_fast(int nevts=0, TString prefix="signal")
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
double GetProb() const
Definition: RhoFitterBase.h:50
void SetType(const TParticlePDG *pdt, int start=0)
TFile * out
Definition: reco_muo.C:20
Double_t M() const
Int_t Remove(RhoCandidate *)
CandList piminus
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