FairRoot/PandaRoot
Public Types | Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
PndAnalysisTask Class Reference

#include <PndAnalysisTask.h>

Inheritance diagram for PndAnalysisTask:

Public Types

typedef std::map< Int_t, Float_t > mapper
 

Public Member Functions

 PndAnalysisTask ()
 
 ~PndAnalysisTask ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 

Protected Attributes

int evcount
 
RhoNeutralParticleSelectorneutralSel
 
RhoPlusParticleSelectorplusSel
 
RhoMinusParticleSelectorminusSel
 
RhoMassParticleSelectorphiMSel
 
RhoMassParticleSelectorpi0MSel
 
RhoMassParticleSelectordsMSel
 
RhoSimpleKaonSelectorkSel
 
RhoSimplePionSelectorpiSel
 
TH1F * phimass
 
TH1F * pi0mass
 
TH1F * dsmass
 
TH1F * ds0mass
 
TH1F * ppmass
 
TH1F * nmult
 

Private Member Functions

virtual void SetParContainers ()
 
 ClassDef (PndAnalysisTask, 1)
 

Private Attributes

TClonesArray * fChargedArray
 
TClonesArray * fNeutralArray
 

Detailed Description

Definition at line 23 of file PndAnalysisTask.h.

Member Typedef Documentation

typedef std::map<Int_t, Float_t> PndAnalysisTask::mapper

Definition at line 27 of file PndAnalysisTask.h.

Constructor & Destructor Documentation

PndAnalysisTask::PndAnalysisTask ( )

Default constructor

Definition at line 35 of file PndAnalysisTask.cxx.

35  :
36  FairTask("Panda Analysis Task")
37 {
38 }
PndAnalysisTask::~PndAnalysisTask ( )

Destructor

Definition at line 42 of file PndAnalysisTask.cxx.

42 { }

Member Function Documentation

PndAnalysisTask::ClassDef ( PndAnalysisTask  ,
 
)
private
void PndAnalysisTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 132 of file PndAnalysisTask.cxx.

References RhoCandList::Add(), RhoCandList::Cleanup(), RhoCandList::Combine(), ds0mass, dsmass, dsMSel, evcount, fChargedArray, fNeutralArray, RhoCandList::GetLength(), RhoFactory::Instance(), kSel, RhoCandidate::M(), minusSel, neutralSel, RhoCandListIterator::Next(), nmult, phimass, phiMSel, pi0mass, pi0MSel, piSel, plusSel, ppmass, RhoCandList::RemoveClones(), RhoFactory::Reset(), and RhoCandList::Select().

133 {
134 
136 
137  if (!(++evcount%100)) { cout <<"evt "<<evcount<<endl; }
138 
139  // **** create all the particle lists we'll need for rebuilding the decay tree
140  //
141  RhoCandList neutralCands,chargedCands, plusCands,minusCands;
142 
143  RhoCandList kpCands,kmCands,piCands;
144 
145  RhoCandList phiCands,pi0Cands,dsCands,ds0Cands,ppCands;
146 
147 
148  RhoCandidate* tc;
149 
150  // **** loop over all Candidates and add them to the list allCands
151  //
152  chargedCands.Cleanup();
153  neutralCands.Cleanup();
154 
155  for (Int_t i1=0; i1<fChargedArray->GetEntriesFast(); i1++) {
156  tc = (RhoCandidate*)fChargedArray->At(i1);
157  chargedCands.Add(tc);
158  }
159 
160  for (Int_t i1=0; i1<fNeutralArray->GetEntriesFast(); i1++) {
161  tc = (RhoCandidate*)fNeutralArray->At(i1);
162  neutralCands.Add(tc);
163  }
164 
165  //cout <<"c:"<<chargedCands.GetLength()<<" n:"<<neutralCands.GetLength()<<endl;
166 
167 
168  // **** select all the basic lists
169  //
170  nmult->Fill(neutralCands.GetLength());
171 
172  plusCands.Select(chargedCands ,plusSel);
173  minusCands.Select(chargedCands ,minusSel);
174 
175  // **** pid selection
176  //
177  kpCands.Select(plusCands ,kSel);
178  kmCands.Select(minusCands ,kSel);
179  piCands.Select(chargedCands ,piSel);
180 
181  // **** now start combining all composits; inbetween plot masses
182  // before using the mass selectors
183  //
184  phiCands.Combine(kpCands,kmCands);
185 
186  RhoCandListIterator iterPhi(phiCands);
187  while ((tc=iterPhi.Next())) { phimass->Fill(tc->M()); }
188  phiCands.Select(phiMSel);
189 
190  dsCands.Combine(phiCands,piCands);
191 
192  RhoCandListIterator iterDs(dsCands);
193  while ((tc=iterDs.Next())) { dsmass->Fill(tc->M()); }
194  dsCands.Select(dsMSel);
195 
196  pi0Cands.Combine(neutralCands,neutralCands);
197 
198  RhoCandListIterator iterPi0(pi0Cands);
199  while ((tc=iterPi0.Next())) { pi0mass->Fill(tc->M()); }
200  pi0Cands.Select(pi0MSel);
201 
202 
203  ds0Cands.Combine(dsCands,pi0Cands);
204 
205  RhoCandListIterator iterDs0(ds0Cands);
206  while ((tc=iterDs0.Next())) { ds0mass->Fill(tc->M()); }
207 
208  ppCands.Combine(ds0Cands,dsCands);
209  ppCands.Select(neutralSel);
210 
211  // **** since we didn't care about the D_s charge, combinations might appear in two
212  // different ways; RemoveClones removes double candidates based on the same final states
213  ppCands.RemoveClones();
214 
215  RhoCandListIterator iterPp(ppCands);
216  while ((tc=iterPp.Next())) { ppmass->Fill(tc->M()); }
217 
218 
219 }
TClonesArray * fNeutralArray
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
Int_t RemoveClones()
void Cleanup()
Definition: RhoCandList.cxx:62
RhoSimpleKaonSelector * kSel
Int_t GetLength() const
Definition: RhoCandList.h:46
static void Reset()
Definition: RhoFactory.cxx:28
void Combine(RhoCandList &l1, RhoCandList &l2)
void Select(RhoParticleSelectorBase *pidmgr)
static RhoFactory * Instance()
Definition: RhoFactory.cxx:34
Double_t M() const
RhoMassParticleSelector * phiMSel
RhoMinusParticleSelector * minusSel
TClonesArray * fChargedArray
RhoNeutralParticleSelector * neutralSel
RhoSimplePionSelector * piSel
RhoMassParticleSelector * dsMSel
RhoPlusParticleSelector * plusSel
RhoMassParticleSelector * pi0MSel
void PndAnalysisTask::Finish ( )
virtual

Definition at line 222 of file PndAnalysisTask.cxx.

References ds0mass, dsmass, nmult, phimass, pi0mass, and ppmass.

223 {
224  phimass->Write();
225  pi0mass->Write();
226  dsmass->Write();
227  ds0mass->Write();
228  ppmass->Write();
229  nmult->Write();
230 
231 }
InitStatus PndAnalysisTask::Init ( )
virtual

Virtual method Init

Definition at line 48 of file PndAnalysisTask.cxx.

References ds0mass, dsmass, dsMSel, evcount, fChargedArray, fNeutralArray, kSel, minusSel, neutralSel, nmult, phimass, phiMSel, pi0mass, pi0MSel, piSel, plusSel, ppmass, and RhoParticleSelectorBase::SetCriterion().

49 {
50 
51  //cout << " Inside the Init function****" << endl;
52 
53  //FairDetector::Initialize();
54  //FairRun* sim = FairRun::Instance();
55  //FairRuntimeDb* rtdb=sim->GetRuntimeDb();
56 
57  // Get RootManager
58  FairRootManager* ioman = FairRootManager::Instance();
59  if ( ! ioman ) {
60  cout << "-E- PndEmcHitProducer::Init: "
61  << "RootManager not instantiated!" << endl;
62  return kFATAL;
63  }
64 
65 // Get input array
66  fChargedArray = (TClonesArray*) ioman->GetObject("PndChargedCandidates");
67  fNeutralArray = (TClonesArray*) ioman->GetObject("PndNeutralCandidates");
68 
69  if ( !fChargedArray && !fNeutralArray) {
70  cout << "-W- PndAnalysisTask::Init: "
71  << "No PndChargedCandidates && PndNeutralCandidates array!" << endl;
72  return kERROR;
73  }
74 
75  // Create and register output array
76  cout << "-I- PndAnalysisTask: Intialization successfull" << endl;
77 
78  phimass = new TH1F("phimass","phi cands",100,0.95,1.1);
79  pi0mass = new TH1F("pi0mass","pi0 cands",100,0.135-0.03,0.135+0.03);
80  dsmass = new TH1F("dsmass","Ds cands",100,1.968-0.03,1.968+0.03);
81  ds0mass = new TH1F("ds0mass","Ds0 cands",100,2.317-0.05,2.317+0.05);
82  ppmass = new TH1F("ppmass","pbarp cands",100,4.306-0.1,4.306+0.1);
83 
84  nmult=new TH1F("nmult","# neutrals",15,0,15);
85 
86  phimass->SetMinimum(0);
87  pi0mass->SetMinimum(0);
88  dsmass->SetMinimum(0);
89  ds0mass->SetMinimum(0);
90  ppmass->SetMinimum(0);
91 
92  // **** create and configure the selectors/filters we'd like to use later
93  //
94  //chargedSel = new RhoChargedParticleSelector;
98 
99  // **** mass selectors for the resonances/composites
100  //
101  phiMSel = new RhoMassParticleSelector("phiSelector" , 1.0195 , 0.01);
102  pi0MSel = new RhoMassParticleSelector("pi0Selector" , 0.135 , 0.005);
103  dsMSel = new RhoMassParticleSelector("dsSelector" , 1.9685 , 0.01);
104 
105  kSel = new RhoSimpleKaonSelector();
106  kSel->SetCriterion("loose");
108  piSel->SetCriterion("veryLoose");
109 
110  evcount=0;
111 
112  return kSUCCESS;
113 
114 }
TClonesArray * fNeutralArray
RhoSimpleKaonSelector * kSel
RhoMassParticleSelector * phiMSel
virtual void SetCriterion(const char *crit)
RhoMinusParticleSelector * minusSel
TClonesArray * fChargedArray
RhoNeutralParticleSelector * neutralSel
RhoSimplePionSelector * piSel
RhoMassParticleSelector * dsMSel
RhoPlusParticleSelector * plusSel
RhoMassParticleSelector * pi0MSel
void PndAnalysisTask::SetParContainers ( )
privatevirtual

Geo file to use Get parameter containers

Definition at line 116 of file PndAnalysisTask.cxx.

References run.

117 {
118 
119  // Get run and runtime database
120  FairRun* run = FairRun::Instance();
121  if ( ! run ) { Fatal("SetParContainers", "No analysis run"); }
122 
123  //FairRuntimeDb* db = run->GetRuntimeDb();
124  //if ( ! db ) Fatal("SetParContainers", "No runtime database");
125 
126 
127 }
Int_t run
Definition: autocutx.C:47

Member Data Documentation

TH1F* PndAnalysisTask::ds0mass
protected

Definition at line 72 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().

TH1F* PndAnalysisTask::dsmass
protected

Definition at line 71 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().

RhoMassParticleSelector* PndAnalysisTask::dsMSel
protected

Definition at line 64 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

int PndAnalysisTask::evcount
protected

Definition at line 49 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndAnalysisTask::fChargedArray
private

Input array of TpcLheTrack

Definition at line 80 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndAnalysisTask::fNeutralArray
private

Definition at line 81 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

RhoSimpleKaonSelector* PndAnalysisTask::kSel
protected

Definition at line 65 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

RhoMinusParticleSelector* PndAnalysisTask::minusSel
protected

Definition at line 58 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

RhoNeutralParticleSelector* PndAnalysisTask::neutralSel
protected

Definition at line 56 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

TH1F* PndAnalysisTask::nmult
protected

Definition at line 75 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().

TH1F* PndAnalysisTask::phimass
protected

book all the histograms

Definition at line 69 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().

RhoMassParticleSelector* PndAnalysisTask::phiMSel
protected

Definition at line 62 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

TH1F* PndAnalysisTask::pi0mass
protected

Definition at line 70 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().

RhoMassParticleSelector* PndAnalysisTask::pi0MSel
protected

Definition at line 63 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

RhoSimplePionSelector* PndAnalysisTask::piSel
protected

Definition at line 66 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

RhoPlusParticleSelector* PndAnalysisTask::plusSel
protected

Definition at line 57 of file PndAnalysisTask.h.

Referenced by Exec(), and Init().

TH1F* PndAnalysisTask::ppmass
protected

Definition at line 73 of file PndAnalysisTask.h.

Referenced by Exec(), Finish(), and Init().


The documentation for this class was generated from the following files: