FairRoot/PandaRoot
Functions
rho/tut_ana_comb.C File Reference

Go to the source code of this file.

Functions

void tut_ana_comb (int nevts=0, TString prefix="signal")
 

Function Documentation

void tut_ana_comb ( int  nevts = 0,
TString  prefix = "signal" 
)

Definition at line 8 of file rho/tut_ana_comb.C.

References RhoCandList::Append(), RhoCandList::Cleanup(), RhoCandList::Combine(), PndAnalysis::FillList(), fRun, PndAnalysis::GetEntries(), PndAnalysis::GetEvent(), RhoCandList::GetLength(), i, muminus, muplus, out, P, piminus, piplus, PndAnalysis::PndAnalysis(), rtdb, RhoCandList::Select(), and TString.

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  fRun->SetUseFairLinks(kTRUE);
30 
31  // *** setup parameter database
32  FairParRootFileIo* parIO = new FairParRootFileIo();
33  parIO->open(inParFile);
34  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
35  parIOPid->open(pidParFile.Data(),"in");
36 
37  rtdb->setFirstInput(parIO);
38  rtdb->setSecondInput(parIOPid);
39  rtdb->setOutput(parIO);
40 
41  fRun->SetOutputFile(OutFile);
42  fRun->Init();
43 
44  // *** create an output file for all histograms
45  TFile *out = TFile::Open(prefix+"_ana_comb.root","RECREATE");
46 
47  // *** create some histograms
48  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
49  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
50 
51  TH1F *hjpsim_pcut = new TH1F("hjpsim_pcut","J/#psi mass (comb. by hand with p cut)",200,0,4.5);
52  TH1F *hpsim_pcut = new TH1F("hpsim_pcut","#psi(2S) mass (comb. by hand with p cut))",200,0,5);
53 
54  //
55  // Now the analysis stuff comes...
56  //
57 
58 
59  // *** the data reader object
60  PndAnalysis* theAnalysis = new PndAnalysis();
61  if (nevts==0) nevts= theAnalysis->GetEntries();
62 
63  // *** RhoCandLists for the analysis
64  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
65 
66  // *** Mass selector for the jpsi cands
67  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",3.096,1.0);
68 
69  // ***
70  // the event loop
71  // ***
72  while (theAnalysis->GetEvent() && i++<nevts)
73  {
74  if ((i%100)==0) cout<<"evt " << i << endl;
75 
76  // *** Select with no PID info ('All'); type and mass are set
77  theAnalysis->FillList(muplus, "MuonAllPlus");
78  theAnalysis->FillList(muminus, "MuonAllMinus");
79  theAnalysis->FillList(piplus, "PionAllPlus");
80  theAnalysis->FillList(piminus, "PionAllMinus");
81 
82  std::cout << "MuPlus:" << std::endl;
83  for (int k = 0; k < muplus.GetLength(); k++){
84  std::cout << k << " : ";
85  ((FairMultiLinkedData_Interface*)muplus[k])->Print();
86  std::cout << std::endl;
87  }
88 
89  std::cout << "MuMinus:" << std::endl;
90  for (int k = 0; k < muminus.GetLength(); k++){
91  std::cout << k << " : ";
92  ((FairMultiLinkedData_Interface*)muminus[k])->Print();
93  std::cout << std::endl;
94  }
95 
96  // ***
97  // *** SIMPLE COMBINATORICS for J/psi -> mu+ mu-
98  // ***
99  jpsi.Combine(muplus, muminus);
100  for (j=0;j<jpsi.GetLength();++j) hjpsim_all->Fill( jpsi[j]->M() );
101 
102  std::cout << "JPsi:" << std::endl;
103  for (int k = 0; k < jpsi.GetLength(); k++){
104  std::cout << k << " : ";
105  ((FairMultiLinkedData_Interface*)jpsi[k])->Print();
106  std::cout << std::endl;
107  }
108 
109 
110  // *** some rough mass selection
111  jpsi.Select(jpsiMassSel);
112 
113  // ***
114  // *** SIMPLE COMBINATORICS for psi(2S) -> J/psi pi+ pi-
115  // ***
116  psi2s.Combine(jpsi, piplus, piminus);
117  for (j=0;j<psi2s.GetLength();++j) hpsim_all->Fill( psi2s[j]->M() );
118 
119  // ***
120  // *** COMBINATORICS BY HAND for J/psi -> mu+ mu- with extra p cut
121  // ***
122 
123  // clean jpsi list and fill it with mu+ mu- combinations
124  jpsi.Cleanup();
125 
126  for (j=0;j<muplus.GetLength();++j)
127  {
128  if (muplus[j]->P()<0.3) continue; // apply momentum selection (can be done easier
129  // with a momentum selector, just for demonstration here)
130  for (k=0;k<muminus.GetLength();++k)
131  {
132  if (muminus[k]->P()<0.3) continue;
133 
134  RhoCandidate *combCand = muplus[j]->Combine(muminus[k]);
135  jpsi.Append(combCand);
136  }
137  }
138 
139  for (j=0;j<jpsi.GetLength();++j) hjpsim_pcut->Fill( jpsi[j]->M() );
140 
141  // *** some rough mass selection
142  jpsi.Select(jpsiMassSel);
143 
144  // ***
145  // *** COMBINATORICS BY HAND for psi(2S) -> J/psi pi+ pi- with extra p cut
146  // ***
147 
148  // clean psi(2S) list and fill it with J/psi pi+ mpi- combinations
149  psi2s.Cleanup();
150 
151  for (j=0;j<jpsi.GetLength();++j)
152  {
153  for (k=0;k<piplus.GetLength();++k)
154  {
155  // momentum selection and check clash with J/psi
156  if ( piplus[k]->P()<0.2 || piplus[k]->Overlaps(jpsi[j]) ) continue;
157 
158  for (l=0;l<piminus.GetLength();++l)
159  {
160  // momentum selection and check clash with J/psi
161  if ( piminus[l]->P()<0.2 || piminus[l]->Overlaps(jpsi[j]) ) continue;
162 
163  RhoCandidate *combCand = jpsi[j]->Combine(piplus[k], piminus[l]);
164  psi2s.Append(combCand);
165  }
166  }
167  }
168 
169  for (j=0;j<psi2s.GetLength();++j) hpsim_pcut->Fill( psi2s[j]->M() );
170 
171 
172 
173  }
174 
175  // *** write out all the histos
176  out->cd();
177 
178  hjpsim_all->Write();
179  hpsim_all->Write();
180 
181  hjpsim_pcut->Write();
182  hpsim_pcut->Write();
183 
184  out->Save();
185 
186 }
Int_t GetEntries()
void Append(const RhoCandidate *c)
Definition: RhoCandList.h:52
void Cleanup()
Definition: RhoCandList.cxx:62
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)
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
TFile * out
Definition: reco_muo.C:20
GeV c P
CandList piminus
Int_t GetEvent(Int_t n=-1)
CandList muminus