FairRoot/PandaRoot
tut_ana_d0_qahelper.C
Go to the documentation of this file.
1 class FairRunAna;
2 
3 // Check if file exists and readable
4 bool checkfile(TString fn)
5 {
6  bool fileok=true;
7  TFile fff(fn, "READ");
8  if (fff.IsZombie()) fileok=false;
9  TTree *t=(TTree*)fff.Get("pndsim");
10  if (t==0x0) fileok=false;
11 
12  if (!fileok) cout <<"Skipping broken file '"<<fn<<"'"<<endl;
13  fff.Close();
14 
15  return fileok;
16 }
17 
18 // Add many input files to FairRunAna
19 void attachFiles(FairRunAna* fRun, TString pref, int min, int max)
20 {
21  bool firstfile=true;
22  for (int i=min;i<=max;++i) {
23  TString fname = TString::Format("%s_%d_pid.root",pref.Data(),i);
24 
25  if (checkfile(fname) ) {
26  if (firstfile) fRun->SetInputFile(fname);
27  else fRun->AddFile(fname);
28  firstfile=false;
29  }
30  }
31 }
32 
33 void tut_ana_d0_qahelper(TString pref="pid_complete.root", int min=-1, int max=1, int nevts=0)
34 {
35  TString OutFile, inParFile;
36 
37  double pbarmom = 9.808065;
38  double mp=0.938272;
39  TLorentzVector ini(0,0,pbarmom, sqrt(pbarmom*pbarmom+mp*mp)+mp);
40 
41  // *** add input files
42  FairRunAna *fRun= new FairRunAna();
43 
44  // *** just open one file
45  if (min<0)
46  {
47  // *** set input file, output file and par file
48  OutFile = "ana.root";
49  inParFile = "simparams.root";
50 
51  fRun->SetInputFile(pref);
52  }
53  // *** attach many input files
54  else {
55  // *** set output file and par file
56  OutFile = TString::Format("%s_ana_%d_%d.root",pref.Data(), min, max);
57  inParFile = TString::Format("%s_%d_par.root",pref.Data(), min);
58 
59  attachFiles(fRun, pref, min, max);
60  }
61 
62  // *** PID table with selection thresholds; can be modified by the user
63  TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
64 
65  // *** initialization
66  FairLogger::GetLogger()->SetLogToFile(kFALSE);
67  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
68 
69  // *** setup parameter database
70  FairParRootFileIo* parIO = new FairParRootFileIo();
71  parIO->open(inParFile);
72  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
73  parIOPid->open(pidParFile.Data(),"in");
74 
75  rtdb->setFirstInput(parIO);
76  rtdb->setSecondInput(parIOPid);
77  rtdb->setOutput(parIO);
78 
79  fRun->SetOutputFile("dummyout.root");
80  fRun->Init();
81 
82  // *** create an output file for all histograms
83  TFile *out = TFile::Open(OutFile,"RECREATE");
84 
85  // *** the data reader object
86  PndAnalysis* theAnalysis = new PndAnalysis();
87  if (nevts==0) nevts= theAnalysis->GetEntries();
88 
89  int i=0;
90 
91  RhoCandList kp, km, pip, pim, d0, d0b, all;
92  TString pidsel = "PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts";
93 
94  RhoMassParticleSelector *fD0Sel = new RhoMassParticleSelector("D0Sel", 1.864, 1.0);
95  RhoVtxPoca fVtxPoca;
96 
97  RhoTuple *nd0 = new RhoTuple("nd0","my D0 tuple");
98 
99  // *** practical helper class to dump stuff to ntuple
100  PndRhoTupleQA qa(theAnalysis,pbarmom);
101 
102  // *** the event loop
103  while (theAnalysis->GetEvent() && i++<nevts)
104  {
105  if ((i%100)==0) cout<<"evt " << i << endl;
106 
107  // ***
108  // *** HERE YOUR ANALYSIS CODE GOES!
109  // ***
110  theAnalysis->FillList( all, "All", pidsel); // only needed for PndEventShape object
111  theAnalysis->FillList( kp, "KaonAllPlus", pidsel);
112  theAnalysis->FillList( km, "KaonAllMinus", pidsel);
113  theAnalysis->FillList( pip, "PionAllPlus", pidsel);
114  theAnalysis->FillList( pim, "PionAllMinus", pidsel);
115 
116  PndEventShape evsh(all,ini,0.03,0.1);
117 
118  d0.Combine(km, pip);
119  d0.SetType(421);
120  d0b.Combine(kp, pim);
121  d0b.SetType(-421);
122  d0.Append(d0b);
123 
124  // *** do some preselection
125  d0.Select(fD0Sel);
126 
127  // *** store all kinf of information to RhoTuple
128  for (int j=0;j<d0.GetLength();++j)
129  {
130  nd0->Column("ev", (Float_t) i);
131  nd0->Column("cand", (Float_t) j);
132  nd0->Column("d0miss", (Float_t) (ini-(d0[j]->P4())).M());
133 
134  qa.qaComp("d0", d0[j], nd0);
135 
136  qa.qaP4("beam", ini, nd0);
137  qa.qaEventShapeShort("es",&evsh,nd0);
138 
139  nd0->DumpData();
140  }
141 
142 
143  }
144 
145  // *** write out all the histos
146  out->cd();
147  nd0->GetInternalTree()->Write();
148  out->Save();
149  out->Close();
150 
151 }
152 
153 
154 
155 
156 
157 
158 
159 
160 
161 
162 
163 
164 
165 
Int_t GetEntries()
Double_t x0
Definition: checkhelixhit.C:70
void Append(const RhoCandidate *c)
Definition: RhoCandList.h:52
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
void qaP4(TString pre, TLorentzVector c, RhoTuple *n, bool skip=false)
static const double mp
Definition: mzparameters.h:11
FairRunAna * fRun
Definition: hit_dirc.C:58
Double_t d0
Definition: checkhelixhit.C:59
void Combine(RhoCandList &l1, RhoCandList &l2)
void attachFiles(FairRunAna *fRun, TString pref, int min, int max)
void qaComp(TString pre, RhoCandidate *c, RhoTuple *n, bool covs=false, bool pulls=false)
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
void Column(const char *label, Bool_t value, Bool_t defval=0, const char *block=0)
Definition: RhoTuple.cxx:56
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
void SetType(const TParticlePDG *pdt, int start=0)
TFile * out
Definition: reco_muo.C:20
void DumpData()
Definition: RhoTuple.cxx:391
bool checkfile(TString fn)
void qaEventShapeShort(TString pre, PndEventShape *evsh, RhoTuple *n)
TTree * GetInternalTree()
Definition: RhoTuple.h:207
TTree * t
Definition: bump_analys.C:13
Int_t GetEvent(Int_t n=-1)
void tut_ana_d0_qahelper(TString pref="pid_complete.root", int min=-1, int max=1, int nevts=0)