FairRoot/PandaRoot
Functions
ana_fast_psi2s.C File Reference
#include "TFile.h"
#include "TH1F.h"
#include "TH2F.h"
#include "TCanvas.h"
#include "TStopwatch.h"
#include "TROOT.h"
#include "TSystem.h"
#include "TTree.h"
#include "TString.h"

Go to the source code of this file.

Functions

void printCand (TLorentzVector l, TVector3 p)
 
int ana_fast_psi2s (TString fname="sim_fast.root", int nevts=0)
 

Function Documentation

int ana_fast_psi2s ( TString  fname = "sim_fast.root",
int  nevts = 0 
)

Definition at line 18 of file ana_fast_psi2s.C.

References c1, ctime, Double_t, Rho4CFitter::FitConserveMasses(), i, printf(), rtime, and timer.

19 {
20  TStopwatch timer;
21  timer.Start();
22 
23  TCanvas *c1=new TCanvas("c1","c1",600,600);
24  c1->Divide(2,2);
25 
26  PndEventReader evr(fname);
27  if (nevts==0) nevts=evr.GetEntries();
28 
29  // **** create and setup some histos for QA plots
30  //
31  TH1F *jpsimass = new TH1F("jpsimass","J/psi cands",100,3.1-0.3,3.1+0.3);
32  TH1F *jpsi2mass = new TH1F("jpsi2mass","J/psi cands 4C fit",100,3.1-0.3,3.1+0.3);
33 
34  TH1F *ppmass = new TH1F("ppmass","pbarp cands",100,3.68598-0.3,3.68598+0.3);
35  TH1F *pp2mass = new TH1F("pp2mass","pbarp fitted",100,3.68598-0.02,3.68598+0.02);
36 
37  // **** create all the particle lists we'll need for rebuilding the decay tree
38  //
39  TCandList pip, pim, ep, em;
40  TCandList jpsi,pp;
41 
42  TPidMassSelector *jpsiMSel = new TPidMassSelector("jpsiSelector" , 3.096 , 0.3);
43 
44  TLorentzVector ini(0,0,6.23164,7.24015);
45 
46  int i=0,j=0;
47 
48  // **** loop over all _events_
49  //
50  while (evr.GetEvent() && ++i<nevts)
51  {
52  if (!(i%100))
53  cout <<"evt "<<i<<endl;
54 
55  //cout<<"e mass:" <<TRho::Instance()->GetPDG()->GetParticle(11)->Mass()<<endl;
56 
57  evr.FillList(pip,"PionVeryLoosePlus");
58  evr.FillList(pim,"PionVeryLooseMinus");
59  evr.FillList(ep, "ElectronVeryLoosePlus");
60  evr.FillList(em, "ElectronVeryLooseMinus");
61 
62 /* for (j=0;j<ep.GetLength();++j)
63  {
64  ep[j].SetMass(TRho::Instance()->GetPDG()->GetParticle(11)->Mass());
65  cout <<"e: "<<ep[j].Uid()<<" ";
66 
67  }
68  cout <<endl;
69  for (j=0;j<em.GetLength();++j) cout <<"e-:"<<em[j].Uid()<<" ";
70  cout <<endl;*/
71 
72  jpsi.Combine(ep,em);
73  //jpsi.Select(jpsiMSel);
74  for (j=0;j<jpsi.GetLength();++j)
75  {
76  jpsimass->Fill(jpsi[j].M());
77  }
78 
79  // for the 4C fit we need to add the daughters individually
80  pp.Combine(jpsi,pip,pim);
81 
82  for (j=0;j<pp.GetLength();++j)
83  {
84  ppmass->Fill(pp[j].M());
85 
86  //do the 4C fit on the pbar p System
87  Rho4CFitter fitter(pp[j],ini);
88 
89  fitter.FitConserveMasses();
90 
91  TCandidate *ppfit=fitter.FittedCand(pp[j]);
92  pp2mass->Fill(ppfit->M());
93 
94  TCandidate *epfit=fitter.FittedCand(*(pp[j].Daughter(0)->Daughter(0)));
95  TCandidate *emfit=fitter.FittedCand(*(pp[j].Daughter(0)->Daughter(1)));
96 
97  TLorentzVector sum=epfit->P4()+emfit->P4();
98 
99  jpsi2mass->Fill(sum.Mag());
100 
101  }
102  }
103 
104  // **** plot all that stuff
105  //
106 
107  c1->cd(1); jpsimass->Draw();
108  c1->cd(2); jpsi2mass->Draw();
109  c1->cd(3); ppmass->Draw();
110  c1->cd(4); pp2mass->Draw();
111  c1->cd();
112 
113  timer.Stop();
114  Double_t rtime = timer.RealTime();
115  Double_t ctime = timer.CpuTime();
116 
117  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
118 
119  return 0;
120 }
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
Int_t i
Definition: run_full.C:25
Double_t
TStopwatch timer
Definition: hit_dirc.C:51
c1
Definition: plot_dirc.C:35
Double_t ctime
Definition: hit_dirc.C:114
Double_t rtime
Definition: hit_dirc.C:113
void printCand ( TLorentzVector  l,
TVector3  p 
)

Definition at line 12 of file ana_fast_psi2s.C.

13 {
14  cout <<"vtx("<<p.X()<<"/"<<p.Y()<<"/"<<p.Z()<<") p4("<<l.Px()<<"/"<<l.Py()<<"/"<<l.Pz()<<")="<<l.Vect().Mag()<<endl;
15 }
Double_t p
Definition: anasim.C:58