FairRoot/PandaRoot
Functions
tut_ana_ntp_noqa.C File Reference
#include "auxtut.C"

Go to the source code of this file.

Functions

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

Function Documentation

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

Definition at line 17 of file tut_ana_ntp_noqa.C.

References RhoKinFitter::Add4MomConstraint(), RhoKinFitter::AddMassConstraint(), RhoTuple::Column(), RhoCandList::Combine(), RhoCandidate::Daughter(), RhoTuple::DumpData(), f, PndAnalysis::FillList(), RhoKinFitter::Fit(), RhoFitterBase::Fit(), fRun, RhoFitterBase::GetChi2(), PndAnalysis::GetEntries(), PndAnalysis::GetEvent(), RhoTuple::GetInternalTree(), RhoCandList::GetLength(), RhoCandidate::GetPidInfo(), RhoFitterBase::GetProb(), RhoCandidate::GetRecoCandidate(), i, initrun(), RhoCandidate::M(), mct, PndAnalysis::McTruthMatch(), muminus, mup, muplus, out, P, RhoCandidate::P(), RhoCandidate::P3(), piminus, piplus, RhoCandidate::Pos(), RemoveGeoManager(), RhoTuple::RhoTuple(), RhoCandList::SetType(), and TString.

18 {
19  // *** some variables
20  int i=0,j=0, k=0, l=0;
21  gStyle->SetOptFit(1011);
22 
23  // *** Initialize FairRunAna with defaults
24  TString OutFile="out_dummy.root";
25  FairRunAna* fRun = initrun(prefix, OutFile);
26  fRun->Init();
27 
28  // *** create an output file for all histograms
29  TFile *out = TFile::Open(prefix+"_ana_ntp.root","RECREATE");
30 
31  // *** create ntuples for J/psi and psi(2S)
32  RhoTuple *njpsi = new RhoTuple("njpsi","J/psi Analysis");
33  RhoTuple *npsip = new RhoTuple("npsip","psi' Analysis");
34 
35  // *** create some columns which might not be filled sometimes
36  njpsi->Column("tjpsim", 0.0f, -999.9f);
37  njpsi->Column("tjpsip", 0.0f, -999.9f);
38  njpsi->Column("tjpsitht", 0.0f, -999.9f);
39 
40  npsip->Column("tpsim", 0.0f, -999.9f);
41  npsip->Column("tpsip", 0.0f, -999.9f);
42  npsip->Column("tpsitht", 0.0f, -999.9f);
43 
44  //
45  // Now the analysis stuff comes...
46  //
47 
48 
49  // *** the data reader object
50  PndAnalysis* theAnalysis = new PndAnalysis();
51  if (nevts==0) nevts= theAnalysis->GetEntries();
52 
53  // *** RhoCandLists for the analysis
54  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
55 
56  // *** Mass selector for the jpsi cands
57  double m0_jpsi = TDatabasePDG::Instance()->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
58  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",m0_jpsi,1.0);
59 
60  // *** Pid Selection Algorithms
61  TString pidSelection = "PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts";
62 
63  // *** the lorentz vector of the initial psi(2S)
64  TLorentzVector ini(0, 0, 6.231552, 7.240065);
65 
66  // ***
67  // the event loop
68  // ***
69  while (theAnalysis->GetEvent() && i++<nevts)
70  {
71  if ((i%100)==0) cout<<"evt " << i << endl;
72 
73  // *** Select with no PID info ('All'); type and mass are set
74  theAnalysis->FillList(muplus, "MuonAllPlus", pidSelection);
75  theAnalysis->FillList(muminus, "MuonAllMinus", pidSelection);
76  theAnalysis->FillList(piplus, "PionAllPlus", pidSelection);
77  theAnalysis->FillList(piminus, "PionAllMinus", pidSelection);
78 
79  // *** combinatorics for J/psi -> mu+ mu-
80  jpsi.Combine(muplus, muminus);
81  jpsi.SetType(443);
82 
83  // *** some mass pre selection
84  //jpsi.Select(jpsiMassSel);
85 
86  // ***
87  // *** do all kind of analysis and store in N-tuple
88  // ***
89  for (j=0;j<jpsi.GetLength();++j)
90  {
91  // get daughters
92  RhoCandidate *mup = jpsi[j]->Daughter(0);
93  RhoCandidate *mum = jpsi[j]->Daughter(1);
96 
97  // get truth information
98  bool mct = theAnalysis->McTruthMatch(jpsi[j]);
99  RhoCandidate *true_jpsi = jpsi[j]->GetMcTruth();
100 
101  // perform vertex fitter
102  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
103  vtxfitter.Fit();
104 
105  RhoCandidate *fitvtx_jpsi = jpsi[j]->GetFit();
106  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
107  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
108  TVector3 vtxpos(-999.,-999.,-999.);
109  if (fitvtx_jpsi) vtxpos = fitvtx_jpsi->Daughter(0)->Pos();
110 
111  // perform mass fit
112  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
113  mfitter.AddMassConstraint(m0_jpsi); // add the mass constraint
114  mfitter.Fit(); // do fit
115 
116  RhoCandidate *fitmass_jpsi = jpsi[j]->GetFit();
117  double chi2_mass = mfitter.GetChi2(); // get chi2 of fit
118  double prob_mass = mfitter.GetProb(); // access probability of fit
119 
120  // *** now write ntuple information
121 
122  // *** general event info
123  njpsi->Column("ev", (Float_t) i, -999.9f);
124  njpsi->Column("cand", (Float_t) j, -999.9f);
125 
126  // *** basic J/psi info
127  njpsi->Column("jpsim", (Float_t) jpsi[j]->M(), -999.9f);
128  njpsi->Column("jpsip", (Float_t) jpsi[j]->P(), -999.9f);
129  njpsi->Column("jpsipt", (Float_t) jpsi[j]->P3().Pt(), -999.9f);
130  njpsi->Column("jpsitht", (Float_t) jpsi[j]->P3().Theta(), -999.9f);
131  njpsi->Column("jpsimissm", (Float_t) (ini-(jpsi[j]->P4())).M(), -999.9f);
132 
133  // *** MC truth info
134  njpsi->Column("mct", (Float_t) mct, -999.9f);
135  if (true_jpsi)
136  {
137  njpsi->Column("tjpsim", (Float_t) true_jpsi->M(), -999.9f);
138  njpsi->Column("tjpsip", (Float_t) true_jpsi->M(), -999.9f);
139  njpsi->Column("tjpsitht",(Float_t) true_jpsi->P3().Theta(), -999.9f);
140  }
141 
142  // *** fitting info
143  njpsi->Column("jpsimvtx", (Float_t) fitvtx_jpsi->M(), -999.9f);
144  njpsi->Column("chi2vtx", (Float_t) chi2_vtx, -999.9f);
145  njpsi->Column("probvtx", (Float_t) prob_vtx, -999.9f);
146  njpsi->Column("vtxx", (Float_t) vtxpos.X(), -999.9f);
147  njpsi->Column("vtxy", (Float_t) vtxpos.Y(), -999.9f);
148  njpsi->Column("vtxz", (Float_t) vtxpos.Z(), -999.9f);
149 
150  njpsi->Column("jpsimmass", (Float_t) fitmass_jpsi->M(), -999.9f);
151  njpsi->Column("chi2mass", (Float_t) chi2_mass, -999.9f);
152  njpsi->Column("probmass", (Float_t) prob_mass, -999.9f);
153 
154  // *** kinematic info of daughters
155  njpsi->Column("mupp", (Float_t) mup->P(), -999.9f);
156  njpsi->Column("muppt", (Float_t) mup->P3().Pt(), -999.9f);
157  njpsi->Column("muptht", (Float_t) mup->P3().Theta(), -999.9f);
158 
159  njpsi->Column("mump", (Float_t) mum->P(), -999.9f);
160  njpsi->Column("mumpt", (Float_t) mum->P3().Pt(), -999.9f);
161  njpsi->Column("mumtht", (Float_t) mum->P3().Theta(), -999.9f);
162 
163  // *** PID info of daughters
164  njpsi->Column("muppid", (Float_t) mup->GetPidInfo(1), -999.9f);
165  njpsi->Column("mumpid", (Float_t) mum->GetPidInfo(1), -999.9f);
166 
167  // *** and finally FILL Ntuple
168  njpsi->DumpData();
169 
170  }
171 
172  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
173  psi2s.Combine(jpsi, piplus, piminus);
174  psi2s.SetType(88888);
175 
176  for (j=0;j<psi2s.GetLength();++j)
177  {
178  // get daughters
179  RhoCandidate *jp = psi2s[j]->Daughter(0);
180  RhoCandidate *pip = psi2s[j]->Daughter(1);
181  RhoCandidate *pim = psi2s[j]->Daughter(2);
182 
185 
186  // get truth information
187  bool mct = theAnalysis->McTruthMatch(psi2s[j]);
188  RhoCandidate *true_psi = psi2s[j]->GetMcTruth();
189 
190  // do 4C fit
191  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
192  fitter.Add4MomConstraint(ini); // set 4 constraint
193  fitter.Fit(); // do fit
194  RhoCandidate *fit4c_jpsi = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
195 
196  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
197  double prob_4c = fitter.GetProb(); // access probability of fit
198 
199  // *** general event info
200  npsip->Column("ev", (Float_t) i, -999.9f);
201  npsip->Column("cand", (Float_t) j, -999.9f);
202 
203  // *** basic psi(2s) info
204  npsip->Column("psim", (Float_t) psi2s[j]->M(), -999.9f);
205  npsip->Column("psip", (Float_t) psi2s[j]->P(), -999.9f);
206  npsip->Column("psipt", (Float_t) psi2s[j]->P3().Pt(), -999.9f);
207  npsip->Column("psitht", (Float_t) psi2s[j]->P3().Theta(), -999.9f);
208 
209  // *** basic J/psi info
210  npsip->Column("jpsim", (Float_t) jp->M(), -999.9f);
211  npsip->Column("jpsip", (Float_t) jp->P(), -999.9f);
212  npsip->Column("jpsipt", (Float_t) jp->P3().Pt(), -999.9f);
213  npsip->Column("jpsitht",(Float_t) jp->P3().Theta(), -999.9f);
214 
215  npsip->Column("jpsim4c",(Float_t) fit4c_jpsi->M(), -999.9f);
216 
217  // *** MC truth info
218  npsip->Column("mct", (Float_t) mct, -999.9f);
219  if (true_psi)
220  {
221  npsip->Column("tpsim", (Float_t) true_psi->M(), -999.9f);
222  npsip->Column("tpsip", (Float_t) true_psi->M(), -999.9f);
223  npsip->Column("tpsitht",(Float_t) true_psi->P3().Theta(), -999.9f);
224  }
225 
226  // *** kinematic info of daughters
227  npsip->Column("pipp", (Float_t) pip->P(), -999.9f);
228  npsip->Column("pippt", (Float_t) pip->P3().Pt(), -999.9f);
229  npsip->Column("piptht", (Float_t) pip->P3().Theta(), -999.9f);
230 
231  npsip->Column("pimp", (Float_t) pim->P(), -999.9f);
232  npsip->Column("pimpt", (Float_t) pim->P3().Pt(), -999.9f);
233  npsip->Column("pimtht", (Float_t) pim->P3().Theta(), -999.9f);
234 
235  // *** PID info of daughters
236  npsip->Column("pippid", (Float_t) pip->GetPidInfo(2), -999.9f);
237  npsip->Column("pimpid", (Float_t) pim->GetPidInfo(2), -999.9f);
238 
239  // *** and finally FILL Ntuple
240  npsip->DumpData();
241  }
242 
243  }
244 
245  // *** write out all the histos
246  out->cd();
247 
248  njpsi->GetInternalTree()->Write();
249  npsip->GetInternalTree()->Write();
250 
251  out->Save();
252 
253  // *** temporaty fix to avoid error on macro exit
255 }
Int_t GetEntries()
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
TVector3 Pos() const
Definition: RhoCandidate.h:186
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
PndPidCandidate * GetRecoCandidate() const
Definition: RhoCandidate.h:376
RhoCandidate * Daughter(Int_t n)
CandList piplus
FairRunAna * fRun
Definition: hit_dirc.C:58
CandList muplus
void Combine(RhoCandList &l1, RhoCandList &l2)
TFile * f
Definition: bump_analys.C:12
void Column(const char *label, Bool_t value, Bool_t defval=0, const char *block=0)
Definition: RhoTuple.cxx:56
void SetType(const TParticlePDG *pdt, int start=0)
TFile * out
Definition: reco_muo.C:20
Double_t M() const
Double_t P() const
void DumpData()
Definition: RhoTuple.cxx:391
GeV c P
CandList piminus
TTree * GetInternalTree()
Definition: RhoTuple.h:207
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
Int_t GetEvent(Int_t n=-1)
CandList muminus
PndMCTrack * mct
Definition: hist-t7.C:147
static const double mup
Definition: mzparameters.h:21
double GetPidInfo(int hypo)
TVector3 P3() const
Definition: RhoCandidate.h:199