FairRoot/PandaRoot
rho/tut_ana_mcmatch.C
Go to the documentation of this file.
1 class RhoCandList;
2 class RhoCandidate;
5 class PndAnalysis;
6 
7 
8 void tut_ana_mcmatch(int nevts = 0, TString prefix = "signal")
9 {
10  // *** some variables
11  int i=0,j=0, k=0, l=0;
12  gStyle->SetOptFit(1011);
13 
14  // *** the output file for FairRunAna
15  TString OutFile="out_dummy.root";
16 
17  // *** the files coming from the simulation
18  TString inPidFile = prefix+"_pid.root"; // this file contains the PndPidCandidates and McTruth
19  TString inParFile = prefix+"_par.root";
20 
21  // *** PID table with selection thresholds; can be modified by the user
22  TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
23 
24  // *** initialization
25  FairLogger::GetLogger()->SetLogToFile(kFALSE);
26  FairRunAna* fRun = new FairRunAna();
27  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
28  fRun->SetSource(new FairFileSource(inPidFile));
29 
30  // *** setup parameter database
31  FairParRootFileIo* parIO = new FairParRootFileIo();
32  parIO->open(inParFile);
33  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
34  parIOPid->open(pidParFile.Data(),"in");
35 
36  rtdb->setFirstInput(parIO);
37  rtdb->setSecondInput(parIOPid);
38  rtdb->setOutput(parIO);
39 
40  fRun->SetOutputFile(OutFile);
41  fRun->Init();
42 
43  // *** create an output file for all histograms
44  TFile *out = TFile::Open(prefix+"_ana_mcmatch.root","RECREATE");
45 
46  // *** create some histograms
47  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass",200,0,4.5);
48  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass",200,0,5);
49 
50  TH1F *hjpsim_nm = new TH1F("hjpsim_nm","J/#psi mass (no truth match)",200,0,4.5);
51  TH1F *hpsim_nm = new TH1F("hpsim_nm","#psi(2S) mass (no truth match)",200,0,5);
52 
53  TH1F *hjpsim_ftm = new TH1F("hjpsim_ftm","J/#psi mass (full truth match)",200,0,4.5);
54  TH1F *hpsim_ftm = new TH1F("hpsim_ftm","#psi(2S) mass (full truth match)",200,0,5);
55 
56  TH1F *hjpsim_diff = new TH1F("hjpsim_diff","J/#psi mass diff to truth",100,-2,2);
57  TH1F *hpsim_diff = new TH1F("hpsim_diff","#psi(2S) mass diff to truth",100,-2,2);
58 
59  //
60  // Now the analysis stuff comes...
61  //
62 
63 
64  // *** the data reader object
65  PndAnalysis* theAnalysis = new PndAnalysis();
66  if (nevts==0) nevts= theAnalysis->GetEntries();
67 
68  // *** RhoCandLists for the analysis
69  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
70 
71  // *** Mass selector for the jpsi cands
72  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",3.096,1.0);
73 
74  // ***
75  // the event loop
76  // ***
77  while (theAnalysis->GetEvent() && i++<nevts)
78  {
79  if ((i%100)==0) cout<<"evt " << i << endl;
80 
81  // *** Select with no PID info ('All'); type and mass are set
82  theAnalysis->FillList(muplus, "MuonAllPlus");
83  theAnalysis->FillList(muminus, "MuonAllMinus");
84  theAnalysis->FillList(piplus, "PionAllPlus");
85  theAnalysis->FillList(piminus, "PionAllMinus");
86 
87  // *** combinatorics for J/psi -> mu+ mu-
88  jpsi.Combine(muplus, muminus);
89 
90  // ***
91  // *** do the TRUTH MATCH for jpsi
92  // ***
93  jpsi.SetType(443);
94 
95  for (j=0;j<jpsi.GetLength();++j)
96  {
97  hjpsim_all->Fill( jpsi[j]->M() );
98 
99  if (theAnalysis->McTruthMatch(jpsi[j]))
100  {
101  hjpsim_ftm->Fill( jpsi[j]->M() );
102  hjpsim_diff->Fill( jpsi[j]->GetMcTruth()->M() - jpsi[j]->M() );
103  }
104  else
105  hjpsim_nm->Fill( jpsi[j]->M() );
106  }
107 
108  // *** some rough mass selection
109  jpsi.Select(jpsiMassSel);
110 
111  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
112  psi2s.Combine(jpsi, piplus, piminus);
113 
114  // ***
115  // *** do the TRUTH MATCH for psi(2S)
116  // ***
117  psi2s.SetType(88888);
118 
119  for (j=0;j<psi2s.GetLength();++j)
120  {
121  hpsim_all->Fill( psi2s[j]->M() );
122 
123  if (theAnalysis->McTruthMatch(psi2s[j]))
124  {
125  hpsim_ftm->Fill( psi2s[j]->M() );
126  hpsim_diff->Fill( psi2s[j]->GetMcTruth()->M() - psi2s[j]->M() );
127  }
128  else
129  hpsim_nm->Fill( psi2s[j]->M() );
130  }
131  }
132 
133  // *** write out all the histos
134  out->cd();
135 
136  hjpsim_all->Write();
137  hpsim_all->Write();
138  hjpsim_nm->Write();
139  hpsim_nm->Write();
140  hjpsim_ftm->Write();
141  hpsim_ftm->Write();
142  hjpsim_diff->Write();
143  hpsim_diff->Write();
144 
145  out->Save();
146 
147 }
Int_t GetEntries()
Int_t i
Definition: run_full.C:25
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 Combine(RhoCandList &l1, RhoCandList &l2)
PndAnalysis(TString tname1="", TString tname2="", TString algnamec="PidAlgoIdealCharged", TString algnamen="PidAlgoIdealNeutral")
Definition: PndAnalysis.cxx:48
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
void SetType(const TParticlePDG *pdt, int start=0)
TFile * out
Definition: reco_muo.C:20
CandList piminus
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
Int_t GetEvent(Int_t n=-1)
CandList muminus
void tut_ana_mcmatch(int nevts=0, TString prefix="signal")