FairRoot/PandaRoot
PndAnalysisTaskExample.cxx
Go to the documentation of this file.
1 /******************************************************
2 
3 JPsi->e+e-
4 Reads the TPC tracks and reconstruct the InvariantMass
5 of J/Psi: Dipak
6 *******************************************************/
7 
8 #include "TClonesArray.h"
9 
10 #include "FairRootManager.h"
11 #include "FairRunAna.h"
12 #include "FairRuntimeDb.h"
13 //#include "PndLheCandidate.h"
14 
15 #include "TVector3.h"
16 #include "TH1F.h"
17 
18 #include "FairRun.h"
19 #include "FairRuntimeDb.h"
20 #include "PndAnalysisTaskExample.h"
21 #include <string>
22 #include <iostream>
23 
24 //RHO stuff
25 #include "RhoBase/RhoCandidate.h"
26 #include "RhoBase/RhoCandList.h"
28 #include "RhoBase/RhoFactory.h"
29 
30 
31 using std::cout;
32 using std::endl;
33 
34 
35 // ----- Default constructor -------------------------------------------
37  FairTask("Panda Analysis Task") {
38 }
39 // -------------------------------------------------------------------------
40 
41 // ----- Destructor ----------------------------------------------------
43 // -------------------------------------------------------------------------
44 
45 
46 
47 // ----- Public method Init --------------------------------------------
49 
50  //cout << " Inside the Init function****" << endl;
51 
52  //FairDetector::Initialize();
53  //FairRun* sim = FairRun::Instance();
54  //FairRuntimeDb* rtdb=sim->GetRuntimeDb();
55 
56  // Get RootManager
57  FairRootManager* ioman = FairRootManager::Instance();
58  if ( ! ioman ) {
59  cout << "-E- PndEmcHitProducer::Init: "
60  << "RootManager not instantiated!" << endl;
61  return kFATAL;
62  }
63 
64  // Get input array
65  fChargedArray = (TClonesArray*) ioman->GetObject("PndChargedCandidates");
66  fNeutralArray = (TClonesArray*) ioman->GetObject("PndNeutralCandidates");
67 
68  if ( !fChargedArray && !fNeutralArray) {
69  cout << "-W- PndAnalysisTaskExample::Init: "
70  << "No PndChargedCandidates && PndNeutralCandidates array!" << endl;
71  return kERROR;
72  }
73 
74  // Create and register output array
75  cout << "-I- PndAnalysisTaskExample: Intialization successfull" << endl;
76 
77  phimass = new TH1F("phimass","phi cands",100,0.95,1.1);
78  pi0mass = new TH1F("pi0mass","pi0 cands",100,0.135-0.03,0.135+0.03);
79  dsmass = new TH1F("dsmass","Ds cands",100,1.968-0.03,1.968+0.03);
80  ds0mass = new TH1F("ds0mass","Ds0 cands",100,2.317-0.05,2.317+0.05);
81  ppmass = new TH1F("ppmass","pbarp cands",100,4.306-0.1,4.306+0.1);
82 
83  nmult=new TH1F("nmult","# neutrals",15,0,15);
84 
85  phimass->SetMinimum(0);
86  pi0mass->SetMinimum(0);
87  dsmass->SetMinimum(0);
88  ds0mass->SetMinimum(0);
89  ppmass->SetMinimum(0);
90 
91  // **** create and configure the selectors/filters we'd like to use later
92  //
93  //chargedSel = new RhoChargedParticleSelector;
97 
98  // **** mass selectors for the resonances/composites
99  //
100  phiMSel = new RhoMassParticleSelector("phiSelector" , 1.0195 , 0.01);
101  pi0MSel = new RhoMassParticleSelector("pi0Selector" , 0.135 , 0.005);
102  dsMSel = new RhoMassParticleSelector("dsSelector" , 1.9685 , 0.01);
103 
104  kSel = new RhoSimpleKaonSelector();
105  kSel->SetCriterion("loose");
107  piSel->SetCriterion("veryLoose");
108 
109  evcount=0;
110 
111  return kSUCCESS;
112 
113 }
114 
116 
117  // Get run and runtime database
118  FairRunAna* run = FairRunAna::Instance();
119  if ( ! run ) Fatal("SetParContainers", "No analysis run");
120 
121  //FairRuntimeDb* db = run->GetRuntimeDb();
122  //if ( ! db ) Fatal("SetParContainers", "No runtime database");
123 
124 
125 }
126 
127 // -------------------------------------------------------------------------
128 
129 // ----- Public method Exec --------------------------------------------
131 
133 
134  if (!(++evcount%100)) cout <<"evt "<<evcount<<endl;
135 
136  // **** create all the particle lists we'll need for rebuilding the decay tree
137  //
138  RhoCandList neutralCands,chargedCands, plusCands,minusCands;
139 
140  RhoCandList kpCands,kmCands,piCands;
141 
142  RhoCandList phiCands,pi0Cands,dsCands,ds0Cands,ppCands;
143 
144 
145  RhoCandidate *tc;
146 
147  // **** loop over all Candidates and add them to the list allCands
148  //
149  chargedCands.Cleanup();
150  neutralCands.Cleanup();
151 
152  for (Int_t i1=0; i1<fChargedArray->GetEntriesFast(); i1++){
153  tc = (RhoCandidate *)fChargedArray->At(i1);
154  chargedCands.Add(*tc);
155  }
156 
157  for (Int_t i1=0; i1<fNeutralArray->GetEntriesFast(); i1++){
158  tc = (RhoCandidate *)fNeutralArray->At(i1);
159  neutralCands.Add(*tc);
160  }
161 
162  //cout <<"c:"<<chargedCands.GetLength()<<" n:"<<neutralCands.GetLength()<<endl;
163 
164 
165  // **** select all the basic lists
166  //
167  nmult->Fill(neutralCands.GetLength());
168 
169  plusCands.Select(chargedCands ,plusSel);
170  minusCands.Select(chargedCands ,minusSel);
171 
172  // **** pid selection
173  //
174  kpCands.Select(plusCands ,kSel);
175  kmCands.Select(minusCands ,kSel);
176  piCands.Select(chargedCands ,piSel);
177 
178  // **** now start combining all composits; inbetween plot masses
179  // before using the mass selectors
180  //
181  phiCands.Combine(kpCands,kmCands);
182 
183  RhoCandListIterator iterPhi(phiCands);
184  while (tc=iterPhi.Next()) phimass->Fill(tc->M());
185  phiCands.Select(phiMSel);
186 
187  dsCands.Combine(phiCands,piCands);
188 
189  RhoCandListIterator iterDs(dsCands);
190  while (tc=iterDs.Next()) dsmass->Fill(tc->M());
191  dsCands.Select(dsMSel);
192 
193  pi0Cands.Combine(neutralCands,neutralCands);
194 
195  RhoCandListIterator iterPi0(pi0Cands);
196  while (tc=iterPi0.Next()) pi0mass->Fill(tc->M());
197  pi0Cands.Select(pi0MSel);
198 
199 
200  ds0Cands.Combine(dsCands,pi0Cands);
201 
202  RhoCandListIterator iterDs0(ds0Cands);
203  while (tc=iterDs0.Next()) ds0mass->Fill(tc->M());
204 
205  ppCands.Combine(ds0Cands,dsCands);
206  ppCands.Select(neutralSel);
207 
208  // **** since we didn't care about the D_s charge, combinations might appear in two
209  // different ways; RemoveClones removes double candidates based on the same final states
210  ppCands.RemoveClones();
211 
212  RhoCandListIterator iterPp(ppCands);
213  while (tc=iterPp.Next()) ppmass->Fill(tc->M());
214 
215 
216 }
217 // -------------------------------------------------------------------------
218 
220 {
221  phimass->Write();
222  pi0mass->Write();
223  dsmass->Write();
224  ds0mass->Write();
225  ppmass->Write();
226  nmult->Write();
227 
228 }
229 
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
Int_t RemoveClones()
Int_t run
Definition: autocutx.C:47
void Cleanup()
Definition: RhoCandList.cxx:62
RhoMinusParticleSelector * minusSel
RhoSimpleKaonSelector * kSel
Int_t GetLength() const
Definition: RhoCandList.h:46
static void Reset()
Definition: RhoFactory.cxx:28
RhoNeutralParticleSelector * neutralSel
RhoPlusParticleSelector * plusSel
void Combine(RhoCandList &l1, RhoCandList &l2)
void Select(RhoParticleSelectorBase *pidmgr)
RhoMassParticleSelector * pi0MSel
static RhoFactory * Instance()
Definition: RhoFactory.cxx:34
Double_t M() const
virtual void SetCriterion(const char *crit)
virtual void Exec(Option_t *opt)
RhoSimplePionSelector * piSel
ClassImp(PndAnaContFact)
RhoMassParticleSelector * dsMSel
RhoMassParticleSelector * phiMSel