FairRoot/PandaRoot
Functions
macro/examples/fastsim/ana_dsinc.C File Reference
#include "TCanvas.h"
#include "TTree.h"

Go to the source code of this file.

Functions

int ana_dsinc (TString fname, int nevts=0)
 

Function Documentation

int ana_dsinc ( TString  fname,
int  nevts = 0 
)

Definition at line 4 of file macro/examples/fastsim/ana_dsinc.C.

References basiclibs(), c1, cm, PndEventInfo::CmFrame(), ctime, Double_t, i, phi, pi, printf(), rtime, and timer.

5 {
6  gROOT->Reset();
7  TStopwatch timer;
8  timer.Start();
9 
10  gROOT->LoadMacro("$VMCWORKDIR/gconfig/basiclibs.C");basiclibs();
11  gSystem->Load("libRho");
12 
13  TCanvas *c1=new TCanvas("c1","c1",600,600);
14  c1->Divide(2,2);
15 
16  TH1F *mphi=new TH1F("mphi","m(phi)",100,0.98,1.06);
17  TH1F *mds=new TH1F("meta","m(ds)",100,1.9685-0.10,1.9685+0.10);
18  TH1F *miss=new TH1F("miss","miss mass",200,2.317-0.2,2.317+0.2);
19  TH1F *sum=new TH1F("sum","m(ds)+miss",200,1.9685+2.317-0.05,1.9685+2.317+0.01);
20 
21  PndEventReader evr(fname);
22 
23  if (nevts==0) nevts=evr.GetEntries();
24  int i=0,j;
25 
26  TCandList kp,km,pi,phi,ds;
27 
28  TPidMassSelector *phiMassSel=new TPidMassSelector("phi",1.0195,0.03);
29 
30  while (evr.GetEvent() && i++<nevts)
31  {
32  if (0==i%100) cout <<"evt "<<i<<endl;
33 
34  PndEventInfo *evtinfo=evr.GetEventInfo();
35 
36  TLorentzVector cm=evtinfo->CmFrame();
37 
38  evr.FillList(kp,"KaonVeryLoosePlus");
39  evr.FillList(km,"KaonVeryLooseMinus");
40  evr.FillList(pi,"PionVeryLoose");
41 
42  phi.Combine(kp,km);
43  phi.Select(phiMassSel);
44  for (j=0;j<phi.GetLength();++j) mphi->Fill(phi[j].M());
45 
46  ds.Combine(phi,pi);
47  for (j=0;j<ds.GetLength();++j)
48  {
49  double dsmass=ds[j].M();
50  mds->Fill(dsmass);
51  double missmass=(cm-ds[j].P4()).M();
52  miss->Fill(missmass);
53  sum->Fill(missmass+dsmass);
54  }
55  }
56 
57  c1->cd(1); mphi->Draw();
58  c1->cd(2); mds->Draw();
59  c1->cd(3); miss->Draw();
60  c1->cd(4); sum->Draw();
61  c1->cd();
62 
63  timer.Stop();
64  Double_t rtime = timer.RealTime();
65  Double_t ctime = timer.CpuTime();
66 
67  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
68  return 0;
69 
70 }
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
basiclibs()
Int_t i
Definition: run_full.C:25
Double_t const cm
#define pi
Definition: createSTT.C:60
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
const TLorentzVector & CmFrame() const
Definition: PndEventInfo.h:47