FairRoot/PandaRoot
macro/examples/fastsim/ana_dsdsj.C
Go to the documentation of this file.
1 int ana_dsdsj(TString fname="dsds_10k.evt.root",int num=0)
2 {
3  TStopwatch timer;
4  timer.Start();
5 
6  gROOT->LoadMacro("$VMCWORKDIR/gconfig/basiclibs.C");
7  basiclibs();
8 
9  // **** the library is needed in order to use the TCandidate and TCandList etc
10  gSystem->Load("libRho");
11 
12  TCanvas *c1=new TCanvas("c1","c1",900,600);
13  c1->Divide(3,2);
14 
15  // **** access the datastructure holding the particle lists
16  //
17  TFile* f = new TFile(fname.Data());
18  TTree *t=f->Get("pndsim") ;
19 
20  // **** for every event, a TCLonesArray with the candidates is stored in pndsim
21  //
22  TClonesArray *fCands=new TClonesArray("TCandidate");
23 
24  t->SetBranchAddress("PndCandidates",&fCands);
25 
26  TCandidate *tc;
27 
28  // **** create and setup some histos for QA plots
29  //
30  TH1F *phimass = new TH1F("phimass","phi cands",100,0.95,1.1);
31  TH1F *pi0mass = new TH1F("pi0mass","pi0 cands",100,0.135-0.03,0.135+0.03);
32  TH1F *dsmass = new TH1F("dsmass","Ds cands",100,1.968-0.03,1.968+0.03);
33  TH1F *ds0mass = new TH1F("ds0mass","Ds0 cands",100,2.317-0.05,2.317+0.05);
34  TH1F *ppmass = new TH1F("ppmass","pbarp cands",100,4.306-0.1,4.306+0.1);
35 
36  TH1F *nmult=new TH1F("nmult","# neutrals",15,0,15);
37 
38  phimass->SetMinimum(0);
39  pi0mass->SetMinimum(0);
40  dsmass->SetMinimum(0);
41  ds0mass->SetMinimum(0);
42  ppmass->SetMinimum(0);
43 
44  if (num==0) num= t->GetEntriesFast();
45  cout <<"\n####### Processing "<<num <<" events...\n"<<endl;
46 
47  // **** create all the particle lists we'll need for rebuilding the decay tree
48  //
49  TCandList allCands,neutralCands,chargedCands, plusCands,minusCands;
50 
51  TCandList kpCands,kmCands,piCands;
52 
53  TCandList phiCands,pi0Cands,dsCands,ds0Cands,ppCands;
54 
55  // **** create and configure the selectors/filters we'd like to use later
56  //
57  TPidChargedSelector *chargedSel = new TPidChargedSelector;
58  TPidNeutralSelector *neutralSel = new TPidNeutralSelector;
59  TPidPlusSelector *plusSel = new TPidPlusSelector;
60  TPidMinusSelector *minusSel = new TPidMinusSelector;
61 
62  // **** mass selectors for the resonances/composites
63  //
64  TPidMassSelector *phiMSel = new TPidMassSelector("phiSelector" , 1.0195 , 0.01);
65  TPidMassSelector *pi0MSel = new TPidMassSelector("pi0Selector" , 0.135 , 0.005);
66  TPidMassSelector *dsMSel = new TPidMassSelector("dsSelector" , 1.9685 , 0.01);
67 
68  TPidSimpleKaonSelector *kSel = new TPidSimpleKaonSelector();
69  kSel->SetCriterion("tight");
70  TPidSimplePionSelector *piSel = new TPidSimplePionSelector();
71  piSel->SetCriterion("loose");
72 
73  int i2;
74 
75  // **** loop over all _events_
76  //
77  for (Int_t j=0; j< num;j++){
78 
79  t->GetEntry(j);
80  TFactory::Instance()->Reset();
81  allCands.Cleanup();
82 
83  // **** loop over all Candidates and add them to the list allCands
84  //
85  for (Int_t i1=0; i1<fCands->GetEntriesFast(); i1++){
86  tc = (TCandidate *)fCands->At(i1);
87  allCands.Add(*tc);
88  }
89 
90  // **** select all the basic list
91  //
92  chargedCands.Select(allCands, chargedSel);
93  neutralCands.Select(allCands, neutralSel);
94 
95  nmult->Fill(neutralCands.GetLength());
96 
97  plusCands.Select(chargedCands ,plusSel);
98  minusCands.Select(chargedCands ,minusSel);
99 
100  // **** pid selection
101  //
102  kpCands.Select(plusCands ,kSel);
103  kmCands.Select(minusCands ,kSel);
104  piCands.Select(chargedCands ,piSel);
105 
106  // **** now start combining all composits; inbetween plot masses
107  // before using the mass selectors
108  //
109  phiCands.Combine(kpCands,kmCands);
110 
111  TCandListIterator iterPhi(phiCands);
112  while (tc=iterPhi.Next()) phimass->Fill(tc->M());
113  phiCands.Select(phiMSel);
114 
115  dsCands.Combine(phiCands,piCands);
116 
117  TCandListIterator iterDs(dsCands);
118  while (tc=iterDs.Next()) dsmass->Fill(tc->M());
119  dsCands.Select(dsMSel);
120 
121  pi0Cands.Combine(neutralCands,neutralCands);
122 
123  TCandListIterator iterPi0(pi0Cands);
124  while (tc=iterPi0.Next()) pi0mass->Fill(tc->M());
125  pi0Cands.Select(pi0MSel);
126 
127  ds0Cands.Combine(dsCands,pi0Cands);
128 
129  TCandListIterator iterDs0(ds0Cands);
130  while (tc=iterDs0.Next()) ds0mass->Fill(tc->M());
131 
132  ppCands.Combine(ds0Cands,dsCands);
133  //cout << ds0Cands.GetLength()<<" "<<dsCands.GetLength()<<" "<<ppCands.GetLength()<<endl;
134  ppCands.Select(neutralSel);
135 
136  // **** since we didn't care about the D_s charge, combinations might appear in two
137  // different ways; RemoveClones removes double candidates based on the same final states
138  ppCands.RemoveClones();
139 
140  TCandListIterator iterPp(ppCands);
141  while (tc=iterPp.Next())
142  {
143  ppmass->Fill(tc->M());
144 
145  printRecursive(tc);
146  }
147 
148  int nch=chargedCands.GetLength();
149  int nn=neutralCands.GetLength();
150  int nphi=phiCands.GetLength();
151  int nds=dsCands.GetLength();
152  int nds0=ds0Cands.GetLength();
153  int npp=ppCands.GetLength();
154 
155  if (npp) cout <<"pp :"<<ppCands[0].GetCharge()<<endl;
156 
157  if (!(j%100))
158  {
159  cout <<"evt:"<<j<<endl;
160  if (nch) cout <<"ch :"<<chargedCands[0].GetCharge()<<endl;
161  if (nn) cout <<"ne :"<<neutralCands[0].GetCharge()<<endl;
162  if (nphi) cout <<"phi:"<<phiCands[0].GetCharge()<<endl;
163  if (nds) cout <<"ds :"<<dsCands[0].GetCharge()<<endl;
164  if (nds0) cout <<"ds0:"<<ds0Cands[0].GetCharge()<<endl;
165  if (npp) cout <<"pp :"<<ppCands[0].GetCharge()<<endl;
166  if (nds>0 && nds0>0 &&dsCands[0].GetCharge()!=ds0Cands[0].GetCharge()) cout <<"****"<<endl;
167  cout <<endl;
168  }
169 
170  }
171 
172  // **** plot all that stuff
173  //
174  c1->cd(1);
175  phimass->Draw();
176  c1->cd(2);
177  pi0mass->Draw();
178  c1->cd(3);
179  dsmass->Draw();
180  c1->cd(4);
181  ds0mass->Draw();
182  c1->cd(5);
183  ppmass->Draw();
184  c1->cd(6);
185  nmult->Draw();
186 
187  timer.Stop();
188  Double_t rtime = timer.RealTime();
189  Double_t ctime = timer.CpuTime();
190  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
191  return 0;
192 }
193 
194 /*void printRecursive(TCandidate *tc, Int_t level=0)
195 {
196 
197  int i=0;
198  int nd=tc->NDaughters();
199  if (nd==0) return;
200  cout <<tc->Uid()<<"("<<level<<") -> ";
201  for (i=0;i<nd;i++) cout<<tc->Daughter(i)->Uid()<<" ";
202  cout <<endl;
203  for (i=0;i<nd;i++) printRecursive(tc->Daughter(i),level+1);
204  if (level==0) cout <<endl;
205 
206 }
207 */
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
basiclibs()
int num[96]
Definition: ranlxd.cxx:381
TClonesArray * fCands
Definition: poormantracks.C:27
Double_t
TStopwatch timer
Definition: hit_dirc.C:51
int printRecursive(TCandidate *tc, int level=0)
TFile * f
Definition: bump_analys.C:12
c1
Definition: plot_dirc.C:35
Double_t ctime
Definition: hit_dirc.C:114
TTree * t
Definition: bump_analys.C:13
Double_t rtime
Definition: hit_dirc.C:113
int ana_dsdsj(TString fname="dsds_10k.evt.root", int num=0)