FairRoot/PandaRoot
thailand2017/solution/tut_ana_pid.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_pid(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_pid.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  //
64  // Now the analysis stuff comes...
65  //
66 
67 
68  // *** the data reader object
69  PndAnalysis* theAnalysis = new PndAnalysis();
70  if (nevts==0) nevts= theAnalysis->GetEntries();
71 
72  // *** RhoCandLists for the analysis
73  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
74 
75  // *** Mass selector for the jpsi cands
76  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",3.096,1.0);
77 
78  // ***
79  // the event loop
80  // ***
81  while (theAnalysis->GetEvent() && i++<nevts)
82  {
83  if ((i%100)==0) cout<<"evt " << i << endl;
84 
85  // ***
86  // *** NO PID combinatorics
87  // ***
88 
89  // *** Select with no PID info ('All'); type and mass are set
90  theAnalysis->FillList(muplus, "MuonAllPlus");
91  theAnalysis->FillList(muminus, "MuonAllMinus");
92  theAnalysis->FillList(piplus, "PionAllPlus");
93  theAnalysis->FillList(piminus, "PionAllMinus");
94 
95  // *** combinatorics for J/psi -> mu+ mu-
96  jpsi.Combine(muplus, muminus);
97 
98  for (j=0;j<jpsi.GetLength();++j) hjpsim_all->Fill( jpsi[j]->M() );
99 
100  // *** some rough mass selection
101  jpsi.Select(jpsiMassSel);
102 
103  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
104  psi2s.Combine(jpsi, piplus, piminus);
105  for (j=0;j<psi2s.GetLength();++j) hpsim_all->Fill( psi2s[j]->M() );
106 
107 
108  // ***
109  // *** TRUE PID combinatorics
110  // ***
111 
112  // *** do MC truth match for PID type
113  SelectTruePid(theAnalysis, muplus);
114  SelectTruePid(theAnalysis, muminus);
115  SelectTruePid(theAnalysis, piplus);
116  SelectTruePid(theAnalysis, piminus);
117 
118  // *** all combinatorics again with true PID
119  jpsi.Combine(muplus, muminus);
120  for (j=0;j<jpsi.GetLength();++j) hjpsim_trpid->Fill( jpsi[j]->M() );
121  jpsi.Select(jpsiMassSel);
122 
123  psi2s.Combine(jpsi, piplus, piminus);
124  for (j=0;j<psi2s.GetLength();++j) hpsim_trpid->Fill( psi2s[j]->M() );
125 
126 
127  // ***
128  // *** LOOSE PID combinatorics
129  // ***
130 
131  // *** and again with PidAlgoMvd;PidAlgoStt;PidAlgoDrc and loose selection
132  theAnalysis->FillList(muplus, "MuonLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
133  theAnalysis->FillList(muminus, "MuonLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
134  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
135  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
136 
137  jpsi.Combine(muplus, muminus);
138  for (j=0;j<jpsi.GetLength();++j) hjpsim_lpid->Fill( jpsi[j]->M() );
139  jpsi.Select(jpsiMassSel);
140 
141  psi2s.Combine(jpsi, piplus, piminus);
142  for (j=0;j<psi2s.GetLength();++j) hpsim_lpid->Fill( psi2s[j]->M() );
143 
144 
145  // ***
146  // *** TIGHT PID combinatorics
147  // ***
148 
149  // *** and again with PidAlgoMvd;PidAlgoStt and tight selection
150  theAnalysis->FillList(muplus, "MuonTightPlus", "PidAlgoMdtHardCuts");
151  theAnalysis->FillList(muminus, "MuonTightMinus", "PidAlgoMdtHardCuts");
152  theAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
153  theAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
154 
155  jpsi.Combine(muplus, muminus);
156  for (j=0;j<jpsi.GetLength();++j) hjpsim_tpid->Fill( jpsi[j]->M() );
157  jpsi.Select(jpsiMassSel);
158 
159  psi2s.Combine(jpsi, piplus, piminus);
160  for (j=0;j<psi2s.GetLength();++j) hpsim_tpid->Fill( psi2s[j]->M() );
161 
162  }
163 
164  // *** change to directory where histograms are created
165  out->cd();
166 
167  // *** plot all histos
168  plotmyhistos();
169 
170  // *** write out all the histos to file
171  int nhist = writemyhistos();
172  cout<<"Writing "<<nhist<<" histograms to file"<<endl;
173  out->Save();
174 
175  // *** temporaty fix to avoid error on macro exit
177 }
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
Int_t i
Definition: run_full.C:25
void RemoveGeoManager()
Definition: auxtut.C:11
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
CandList piplus
FairRunAna * fRun
Definition: hit_dirc.C:58
CandList muplus
void tut_ana_pid(int nevts=0, TString prefix="signal")
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)
TFile * out
Definition: reco_muo.C:20
Int_t Remove(RhoCandidate *)
CandList piminus
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
Int_t GetEvent(Int_t n=-1)
CandList muminus
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