FairRoot/PandaRoot
thailand2017/solution/tut_ana.C
Go to the documentation of this file.
1 class RhoCandList;
2 class RhoCandidate;
5 class PndAnalysis;
6 
7 // **** some auxilliary functions in auxtut.C ****
8 // - FairRunAna* initrun(TString prefix, TString outfile, int min=-1, int max=-1) --> Init FairRunAna
9 // - plotmyhistos() --> Plots all histograms in current TDirectory on a autosized canvas
10 // - writemyhistos() --> Writes all histos in current TFile
11 // - fillM(RhoCandList l, TH1* h) --> Fill mass histogram h with masses of candidates in l
12 // - RemoveGeoManager() --> Temporary fix for error on macro exit
13 // **** some auxilliary functions in auxtut.C ****
14 #include "auxtut.C"
15 
16 
17 // *** routine to only keep PID matched candidates in list
19 {
20  int removed = 0;
21 
22  for (int ii=l.GetLength()-1;ii>=0;--ii)
23  {
24  if ( !(ana->McTruthMatch(l[ii])) )
25  {
26  l.Remove(l[ii]);
27  removed++;
28  }
29  }
30 
31  return removed;
32 }
33 
34 
35 void tut_ana(int nevts = 0, TString prefix = "signal")
36 {
37  // *** some variables
38  int i=0,j=0, k=0, l=0;
39  gStyle->SetOptFit(1011);
40 
41  // *** Initialize FairRunAna with defaults
42  TString OutFile="out_dummy.root";
43  FairRunAna* fRun = initrun(prefix, OutFile);
44  fRun->Init();
45 
46  // *** create an output file for all histograms
47  TFile *out = TFile::Open(prefix+"_ana.root","RECREATE");
48 
49  // *** create some histograms
50  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
51  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
52 
53  TH1F *hjpsim_lpid = new TH1F("hjpsim_lpid","J/#psi mass (loose pid)",200,0,4.5);
54  TH1F *hpsim_lpid = new TH1F("hpsim_lpid","#psi(2S) mass (loose pid)",200,0,5);
55 
56  TH1F *hjpsim_tpid = new TH1F("hjpsim_tpid","J/#psi mass (tight pid)",200,0,4.5);
57  TH1F *hpsim_tpid = new TH1F("hpsim_tpid","#psi(2S) mass (tight pid)",200,0,5);
58 
59  TH1F *hjpsim_trpid = new TH1F("hjpsim_trpid","J/#psi mass (true pid)",200,0,4.5);
60  TH1F *hpsim_trpid = new TH1F("hpsim_trpid","#psi(2S) mass (true pid)",200,0,5);
61 
62 
63  TH1F *hjpsim_ftm = new TH1F("hjpsim_ftm","J/#psi mass (full truth match)",200,0,4.5);
64  TH1F *hpsim_ftm = new TH1F("hpsim_ftm","#psi(2S) mass (full truth match)",200,0,5);
65 
66  TH1F *hjpsim_nm = new TH1F("hjpsim_nm","J/#psi mass (no truth match)",200,0,4.5);
67  TH1F *hpsim_nm = new TH1F("hpsim_nm","#psi(2S) mass (no truth match)",200,0,5);
68 
69  TH1F *hjpsim_diff = new TH1F("hjpsim_diff","J/#psi mass diff to truth",100,-2,2);
70  TH1F *hpsim_diff = new TH1F("hpsim_diff","#psi(2S) mass diff to truth",100,-2,2);
71 
72 
73  TH1F *hjpsim_vf = new TH1F("hjpsim_vf","J/#psi mass (vertex fit)",200,0,4.5);
74  TH1F *hjpsim_4cf = new TH1F("hjpsim_4cf","J/#psi mass (4C fit)",200,0,4.5);
75  TH1F *hjpsim_mcf = new TH1F("hjpsim_mcf","J/#psi mass (mass constraint fit)",200,0,4.5);
76 
77  TH1F *hjpsi_chi2_vf = new TH1F("hjpsi_chi2_vf", "J/#psi: #chi^{2} vertex fit",100,0,10);
78  TH1F *hpsi_chi2_4c = new TH1F("hpsi_chi2_4c", "#psi(2S): #chi^{2} 4C fit",100,0,250);
79  TH1F *hjpsi_chi2_mf = new TH1F("hjpsi_chi2_mf", "J/#psi: #chi^{2} mass fit",100,0,10);
80 
81  TH1F *hjpsi_prob_vf = new TH1F("hjpsi_prob_vf", "J/#psi: Prob vertex fit",100,0,1);
82  TH1F *hpsi_prob_4c = new TH1F("hpsi_prob_4c", "#psi(2S): Prob 4C fit",100,0,1);
83  TH1F *hjpsi_prob_mf = new TH1F("hjpsi_prob_mf", "J/#psi: Prob mass fit",100,0,1);
84 
85  TH2F *hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-2,2,100,-2,2);
86 
87  //
88  // Now the analysis stuff comes...
89  //
90 
91 
92  // *** the data reader object
93  PndAnalysis* theAnalysis = new PndAnalysis();
94  if (nevts==0) nevts= theAnalysis->GetEntries();
95 
96  // *** RhoCandLists for the analysis
97  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
98 
99  // *** Mass selector for the jpsi cands
100  double m0_jpsi = TDatabasePDG::Instance()->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
101  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",m0_jpsi,1.0);
102 
103  // *** the lorentz vector of the initial psi(2S)
104  TLorentzVector ini(0, 0, 6.231552, 7.240065);
105 
106  // ***
107  // the event loop
108  // ***
109  while (theAnalysis->GetEvent() && i++<nevts)
110  {
111  if ((i%100)==0) cout<<"evt " << i << endl;
112 
113  // *** Select with no PID info ('All'); type and mass are set
114  theAnalysis->FillList(muplus, "MuonAllPlus");
115  theAnalysis->FillList(muminus, "MuonAllMinus");
116  theAnalysis->FillList(piplus, "PionAllPlus");
117  theAnalysis->FillList(piminus, "PionAllMinus");
118 
119  // *** combinatorics for J/psi -> mu+ mu-
120  jpsi.Combine(muplus, muminus);
121 
122 
123  // ***
124  // *** do the TRUTH MATCH for jpsi
125  // ***
126  jpsi.SetType(443);
127 
128  for (j=0;j<jpsi.GetLength();++j)
129  {
130  hjpsim_all->Fill( jpsi[j]->M() );
131 
132  if (theAnalysis->McTruthMatch(jpsi[j]))
133  {
134  hjpsim_ftm->Fill( jpsi[j]->M() );
135  hjpsim_diff->Fill( jpsi[j]->GetMcTruth()->M() - jpsi[j]->M() );
136  }
137  else
138  hjpsim_nm->Fill( jpsi[j]->M() );
139  }
140 
141  // ***
142  // *** do VERTEX FIT (J/psi)
143  // ***
144  for (j=0;j<jpsi.GetLength();++j)
145  {
146  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
147  vtxfitter.Fit();
148 
149  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
150  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
151  hjpsi_chi2_vf->Fill(chi2_vtx);
152  hjpsi_prob_vf->Fill(prob_vtx);
153 
154  if ( prob_vtx > 0.01 ) // when good enough, fill some histos
155  {
156  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
157  TVector3 jVtx=jfit->Pos(); // and the decay vertex position
158 
159  hjpsim_vf->Fill(jfit->M());
160  hvpos->Fill(jVtx.X(),jVtx.Y());
161  }
162  }
163 
164  // *** some rough mass selection
165  jpsi.Select(jpsiMassSel);
166 
167  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
168  psi2s.Combine(jpsi, piplus, piminus);
169 
170 
171  // ***
172  // *** do the TRUTH MATCH for psi(2S)
173  // ***
174  psi2s.SetType(88888);
175 
176  for (j=0;j<psi2s.GetLength();++j)
177  {
178  hpsim_all->Fill( psi2s[j]->M() );
179 
180  if (theAnalysis->McTruthMatch(psi2s[j]))
181  {
182  hpsim_ftm->Fill( psi2s[j]->M() );
183  hpsim_diff->Fill( psi2s[j]->GetMcTruth()->M() - psi2s[j]->M() );
184  }
185  else
186  hpsim_nm->Fill( psi2s[j]->M() );
187  }
188 
189 
190  // ***
191  // *** do 4C FIT (initial psi(2S) system)
192  // ***
193  for (j=0;j<psi2s.GetLength();++j)
194  {
195  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
196  fitter.Add4MomConstraint(ini); // set 4 constraint
197  fitter.Fit(); // do fit
198 
199  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
200  double prob_4c = fitter.GetProb(); // access probability of fit
201  hpsi_chi2_4c->Fill(chi2_4c);
202  hpsi_prob_4c->Fill(prob_4c);
203 
204  if ( prob_4c > 0.01 ) // when good enough, fill some histo
205  {
206  RhoCandidate *jfit = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
207 
208  hjpsim_4cf->Fill(jfit->M());
209  }
210  }
211 
212 
213  // ***
214  // *** do MASS CONSTRAINT FIT (J/psi)
215  // ***
216  for (j=0;j<jpsi.GetLength();++j)
217  {
218  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
219  mfitter.AddMassConstraint(m0_jpsi); // add the mass constraint
220  mfitter.Fit(); // do fit
221 
222  double chi2_m = mfitter.GetChi2(); // get chi2 of fit
223  double prob_m = mfitter.GetProb(); // access probability of fit
224  hjpsi_chi2_mf->Fill(chi2_m);
225  hjpsi_prob_mf->Fill(prob_m);
226 
227  if ( prob_m > 0.01 ) // when good enough, fill some histo
228  {
229  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
230  hjpsim_mcf->Fill(jfit->M());
231  }
232  }
233 
234 
235  // ***
236  // *** TRUE PID combinatorics
237  // ***
238 
239  // *** do MC truth match for PID type
240  SelectTruePid(theAnalysis, muplus);
241  SelectTruePid(theAnalysis, muminus);
242  SelectTruePid(theAnalysis, piplus);
243  SelectTruePid(theAnalysis, piminus);
244 
245  // *** all combinatorics again with true PID
246  jpsi.Combine(muplus, muminus);
247  fillM(jpsi, hjpsim_trpid);
248  jpsi.Select(jpsiMassSel);
249 
250  psi2s.Combine(jpsi, piplus, piminus);
251  fillM(psi2s, hpsim_trpid);
252 
253 
254  // ***
255  // *** LOOSE PID combinatorics
256  // ***
257 
258  // *** and again with PidAlgoMvd;PidAlgoStt;PidAlgoDrc and loose selection
259  theAnalysis->FillList(muplus, "MuonLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
260  theAnalysis->FillList(muminus, "MuonLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
261  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
262  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
263 
264  jpsi.Combine(muplus, muminus);
265  fillM(jpsi, hjpsim_lpid);
266  jpsi.Select(jpsiMassSel);
267 
268  psi2s.Combine(jpsi, piplus, piminus);
269  fillM(psi2s, hpsim_lpid);
270 
271 
272  // ***
273  // *** TIGHT PID combinatorics
274  // ***
275 
276  // *** and again with PidAlgoMvd;PidAlgoStt and tight selection
277  theAnalysis->FillList(muplus, "MuonTightPlus", "PidAlgoMdtHardCuts");
278  theAnalysis->FillList(muminus, "MuonTightMinus", "PidAlgoMdtHardCuts");
279  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
280  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
281 
282  jpsi.Combine(muplus, muminus);
283  fillM(jpsi, hjpsim_tpid);
284  jpsi.Select(jpsiMassSel);
285 
286  psi2s.Combine(jpsi, piplus, piminus);
287  fillM(psi2s, hpsim_tpid);
288  }
289 
290  // *** change to directory where histograms are created
291  out->cd();
292 
293  // *** plot all histos
294  plotmyhistos();
295 
296  // *** write out all the histos to file
297  int nhist = writemyhistos();
298  cout<<"Writing "<<nhist<<" histograms to file"<<endl;
299  out->Save();
300 
301  // *** temporaty fix to avoid error on macro exit
303 }
Int_t GetEntries()
int SelectTruePid(PndAnalysis *ana, RhoCandList &l)
FairRunAna * initrun(TString prefix, TString outfile, int min=-1, int max=-1)
Definition: QA/auxi.C:32
void AddMassConstraint(double mass)
Int_t i
Definition: run_full.C:25
void RemoveGeoManager()
Definition: auxtut.C:11
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)
void tut_ana(TString pref="pid_complete.root", int min=-1, int max=1, int nevts=0)
Definition: dec13/tut_ana.C:33
CandList piplus
void Add4MomConstraint(TLorentzVector lv)
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 Select(RhoParticleSelectorBase *pidmgr)
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
int fillM(RhoCandList &l, TH1 *h)
Definition: QA/auxi.C:139
int writemyhistos(int maxy=800, double asp=1.1)
Definition: QA/auxi.C:121
void plotmyhistos(std::vector< TH1 * > h, int maxy=700, double asp=1.1)
Definition: QA/auxi.C:62