FairRoot/PandaRoot
Public Member Functions | Private Attributes | List of all members
PndListProvider Class Reference

#include <PndListProvider.h>

Public Member Functions

 PndListProvider (std::string name, std::string pdgType="")
 
 PndListProvider (std::string name, int pdgcode)
 
virtual ~PndListProvider ()
 
std::string GetName ()
 
bool IsGeneric ()
 
bool ToDump ()
 
void SetType (std::string pdgType)
 
void SetType (int pdgcode)
 
void SetGeneric (bool isgeneric=true)
 
void SetToDump (bool todump=true)
 
void SetHisto (TH1F *h)
 
void SetCandList (RhoCandList &cl)
 
void SetMassSelector (double mean, double width)
 
void AddDaughterPointer (PndListProvider *p)
 
void AddDecayProduct (std::string dname)
 
void AddDaughterType (std::string dtype)
 
void AddDaughterType (int pdgcode)
 
int GetType ()
 
int GetDaughterType (int i)
 
std::string GetDecayProdName (int i)
 
int GetNDaughters ()
 
void GetCandList (RhoCandList &tl)
 
int GetNCandidates ()
 
TH1F * GetHisto ()
 
void FillHisto ()
 
void Reset ()
 
void Print ()
 

Private Attributes

std::string fName
 
RhoCandList fOwnList
 
RhoMinusParticleSelectorfMinusSel
 
RhoPlusParticleSelectorfPlusSel
 
RhoMassParticleSelectorfMassSel
 
TH1F * fHisto
 
TParticlePDG * fParticlePDG
 
std::vector< TParticlePDG * > fDaughterPDG
 
ArgVector fDaughterListNames
 
std::vector< PndListProvider * > fDaughterPointers
 
int fNDaughters
 
bool fIsGeneric
 
bool fToDump
 
bool fCandsUptodate
 
bool fHistoFilled
 

Detailed Description

Definition at line 55 of file PndListProvider.h.

Constructor & Destructor Documentation

PndListProvider::PndListProvider ( std::string  name,
std::string  pdgType = "" 
)

Definition at line 51 of file PndListProvider.cxx.

References RhoCandList::Cleanup(), fCandsUptodate, fDaughterListNames, fDaughterPDG, fHisto, fHistoFilled, fIsGeneric, fMassSel, fMinusSel, fName, fNDaughters, fOwnList, fParticlePDG, fPlusSel, fToDump, and name.

52 {
53  fName=name;
54  fNDaughters=0;
55  fIsGeneric=false;
56  fToDump=false;
57  fHistoFilled=false;
58  fCandsUptodate=false;
59 
60  fHisto=0;
61 
62  fDaughterListNames.clear();
63  fDaughterPDG.clear();
64 
65  fOwnList.Cleanup();
66 
69  fMassSel = 0;
70 
71  if (pdgType!="") {
72  fParticlePDG=TDatabasePDG::Instance()->GetParticle(pdgType.c_str());
73  } else {
74  fParticlePDG=0;
75  }
76 }
RhoMinusParticleSelector * fMinusSel
void Cleanup()
Definition: RhoCandList.cxx:62
RhoCandList fOwnList
RhoMassParticleSelector * fMassSel
RhoPlusParticleSelector * fPlusSel
ArgVector fDaughterListNames
std::string fName
TString name
std::vector< TParticlePDG * > fDaughterPDG
TParticlePDG * fParticlePDG
PndListProvider::PndListProvider ( std::string  name,
int  pdgcode 
)

Definition at line 80 of file PndListProvider.cxx.

References RhoCandList::Cleanup(), fCandsUptodate, fDaughterListNames, fDaughterPDG, fHisto, fHistoFilled, fIsGeneric, fMassSel, fMinusSel, fName, fNDaughters, fOwnList, fParticlePDG, fPlusSel, fToDump, and name.

81 {
82  fName=name;
83  fNDaughters=0;
84  fIsGeneric=false;
85  fHistoFilled=false;
86  fToDump=false;
87  fCandsUptodate=false;
88 
89  fHisto=0;
90 
91  fDaughterListNames.clear();
92  fDaughterPDG.clear();
93 
94  fOwnList.Cleanup();
95 
98  fMassSel = 0;
99 
100  if (pdgcode!=0) {
101  fParticlePDG=TDatabasePDG::Instance()->GetParticle(pdgcode);
102  } else {
103  fParticlePDG=0;
104  }
105 }
RhoMinusParticleSelector * fMinusSel
void Cleanup()
Definition: RhoCandList.cxx:62
RhoCandList fOwnList
RhoMassParticleSelector * fMassSel
RhoPlusParticleSelector * fPlusSel
ArgVector fDaughterListNames
std::string fName
TString name
std::vector< TParticlePDG * > fDaughterPDG
TParticlePDG * fParticlePDG
PndListProvider::~PndListProvider ( )
virtual

Definition at line 110 of file PndListProvider.cxx.

References fMassSel, fMinusSel, fPlusSel, and Reset().

111 {
112  delete fMinusSel;
113  delete fPlusSel;
114  delete fMassSel;
115  Reset();
116 }
RhoMinusParticleSelector * fMinusSel
RhoMassParticleSelector * fMassSel
RhoPlusParticleSelector * fPlusSel

Member Function Documentation

void PndListProvider::AddDaughterPointer ( PndListProvider p)

Definition at line 355 of file PndListProvider.cxx.

References fDaughterPointers.

356 {
357  if (p) { fDaughterPointers.push_back(p); }
358 }
std::vector< PndListProvider * > fDaughterPointers
void PndListProvider::AddDaughterType ( std::string  dtype)

Definition at line 134 of file PndListProvider.cxx.

References fDaughterPDG.

135 {
136  fDaughterPDG.push_back(TDatabasePDG::Instance()->GetParticle(dtype.c_str()));
137 }
std::vector< TParticlePDG * > fDaughterPDG
void PndListProvider::AddDaughterType ( int  pdgcode)

Definition at line 141 of file PndListProvider.cxx.

References fDaughterPDG.

142 {
143  fDaughterPDG.push_back(TDatabasePDG::Instance()->GetParticle(pdgcode));
144 }
std::vector< TParticlePDG * > fDaughterPDG
void PndListProvider::AddDecayProduct ( std::string  dname)

Definition at line 158 of file PndListProvider.cxx.

References fDaughterListNames, and fNDaughters.

159 {
160  fDaughterListNames.push_back(std::string(dname));
161  fNDaughters++;
162 }
ArgVector fDaughterListNames
void PndListProvider::FillHisto ( )

Definition at line 371 of file PndListProvider.cxx.

References fHisto, fHistoFilled, fOwnList, RhoCandList::GetLength(), i, and Mass.

372 {
373  if (fHisto && !fHistoFilled) {
374  for (int i=0; i<fOwnList.GetLength(); ++i) { fHisto->Fill(fOwnList[i]->Mass()); }
375  fHistoFilled=true;
376  }
377 }
Int_t i
Definition: run_full.C:25
Int_t GetLength() const
Definition: RhoCandList.h:46
RhoCandList fOwnList
void PndListProvider::GetCandList ( RhoCandList tl)

Definition at line 224 of file PndListProvider.cxx.

References RhoCandidate::AddDaughterLinkSimple(), Bool_t, c, RhoCandList::Cleanup(), RhoCandList::Combine(), Double_t, fCandsUptodate, fDaughterPointers, fHisto, fMassSel, fNDaughters, fOwnList, GetCandList(), RhoCandList::GetLength(), RhoCandList::GetNumberOfTracks(), i, Mass, RhoCandList::Put(), RhoCandList::Select(), and RhoCandidate::SetMarker().

Referenced by GetCandList().

225 {
226  tcl.Cleanup();
227 
228  if (!fCandsUptodate) {
229  // ********* 2 Daughters
230  if (fNDaughters==2) {
231  RhoCandList fD1List,fD2List;
232  fDaughterPointers[0]->GetCandList(fD1List);
234  fDaughterPointers[1]->GetCandList(fD2List);
235  fOwnList.Combine(fD1List,fD2List);
236  } else {
237  fOwnList.Combine(fD1List,fD1List);
238  }
239 
240  //cout <<"d1:"<<fD1List.GetLength()<<" d2:"<<fD1List.GetLength()<<endl;
241  }
242  // ********* 3 Daughters
243  else if (fNDaughters==3) {
244  RhoCandList fDList[3];
245  PndListProvider* lp[3];
246  for (int i=0; i<3; i++) {
247  lp[i]=fDaughterPointers[i];
248  lp[i]->GetCandList(fDList[i]);
249  }
250  //3 different lists
251  if (lp[0]!=lp[1] && lp[1]!=lp[2] && lp[0]!=lp[2]) {
252  RhoCandList tmpL;
253  tmpL.Combine(fDList[0],fDList[1]);
254  fOwnList.Combine(tmpL,fDList[2]);
255  }
256  //at least 2 different lists
257  else if (lp[0]!=lp[1] || lp[1]!=lp[2] || lp[0]!=lp[2]) {
258  int id1=-1,id2=-1;
259  if (lp[0]==lp[1] && lp[0]!=lp[2]) {id1=0; id2=2;}
260  if (lp[0]==lp[2] && lp[0]!=lp[1]) {id1=0; id2=1;}
261  if (lp[1]==lp[2] && lp[0]!=lp[1]) {id1=1; id2=0;}
262 
263  //exact 2 list are the same
264  if (id1!=-1 && id2!=-1) {
265  RhoCandList tmpL;
266  tmpL.Combine(fDList[id1],fDList[id1]);
267  fOwnList.Combine(tmpL,fDList[id2]);
268  }
269  }
270  //all 3 lists are the same
271  else {
272  /*
273  RhoCandList tmpL;
274 
275  tmpL.Combine(fDList[0],fDList[0]);
276  fOwnList.Combine(tmpL,fDList[0]);
277  //fOwnList.RemoveClones();
278  */
279 
280  TLorentzVector vl;
281  Double_t charge;
282  Bool_t nearby = kTRUE;
283  RhoCandList l1;
284  lp[0]->GetCandList(l1);
285 
286  int endpos = l1.GetNumberOfTracks();
287 
288  //combination of a list with itself 3x
289  for (Int_t comb_i=0; comb_i<endpos; ++comb_i) {
290  for (Int_t comb_j=comb_i+1; comb_j<endpos; ++comb_j) {
291  for (Int_t comb_k=comb_j+1; comb_k<endpos; ++comb_k) {
292  if ( !l1[comb_i]->Overlaps( l1[comb_j] )
293  && !l1[comb_j]->Overlaps( l1[comb_k] )
294  && !l1[comb_i]->Overlaps( l1[comb_k] ) ) {
295  vl=l1[comb_i]->P4()+l1[comb_j]->P4()+l1[comb_k]->P4();
296  charge=l1[comb_i]->Charge()+l1[comb_j]->Charge()+l1[comb_k]->Charge();
297  //if (selector) nearby = selector->Accept(l1[comb_i],l2[comb_k]);
298  if (nearby) {
299  //fill list with new candidate
300  RhoCandidate c(vl,charge);
301  c.SetMarker(l1[comb_i]->GetMarker(0)|l1[comb_j]->GetMarker(0)|l1[comb_k]->GetMarker(0),0);
302  c.SetMarker(l1[comb_i]->GetMarker(1)|l1[comb_j]->GetMarker(1)|l1[comb_k]->GetMarker(1),1);
303  c.SetMarker(l1[comb_i]->GetMarker(2)|l1[comb_j]->GetMarker(2)|l1[comb_k]->GetMarker(2),2);
304  c.SetMarker(l1[comb_i]->GetMarker(3)|l1[comb_j]->GetMarker(3)|l1[comb_k]->GetMarker(3),3);
305  //if (selector!=0)
306  //{
307  // c.SetPosition(selector->GetVertex());
308  // c.SetVect(selector->GetMomentum());
309  // c.SetEnergy(c.E());
310  //}
311  // *************** modified by K Goetzen
312  c.AddDaughterLinkSimple(l1[comb_i]);
313  c.AddDaughterLinkSimple(l1[comb_j]);
314  c.AddDaughterLinkSimple(l1[comb_k]);
315  // ****************
316  fOwnList.Put(&c);
317  } //if nearby
318  }// overlap
319  }//for k
320  }//for j
321  }// for i
322 
323  }//else 3x same
324 
325  }
326 
327  if (fHisto)
328  for (int i=0; i<fOwnList.GetLength(); ++i) { fHisto->Fill(fOwnList[i]->Mass()); }
329 
330  if (fMassSel) { fOwnList.Select(fMassSel); }
331 
332  fCandsUptodate=true;
333  }
334 
335 
336  //cout <<fName<<":"<<fParticlePDG->PdgCode()<<endl;
337  //cout <<fOwnList<<endl;
338 
339  tcl=fOwnList;
340 }
Int_t i
Definition: run_full.C:25
Int_t GetLength() const
Definition: RhoCandList.h:46
RhoCandList fOwnList
RhoMassParticleSelector * fMassSel
std::vector< PndListProvider * > fDaughterPointers
void Combine(RhoCandList &l1, RhoCandList &l2)
Double_t
void Select(RhoParticleSelectorBase *pidmgr)
Int_t GetNumberOfTracks() const
Definition: RhoCandList.cxx:72
void Put(const RhoCandidate *, Int_t i=-1)
Definition: RhoCandList.cxx:77
void GetCandList(RhoCandList &tl)
int PndListProvider::GetDaughterType ( int  i)

Definition at line 148 of file PndListProvider.cxx.

References fDaughterPDG, and i.

149 {
150  if ((int)fDaughterPDG.size()>i) {
151  return fDaughterPDG[i]->PdgCode();
152  }
153  return -1;
154 }
Int_t i
Definition: run_full.C:25
std::vector< TParticlePDG * > fDaughterPDG
std::string PndListProvider::GetDecayProdName ( int  i)

Definition at line 346 of file PndListProvider.cxx.

References fDaughterListNames, fNDaughters, and i.

347 {
348  if (i<fNDaughters) {
349  return fDaughterListNames[i];
350  } else { return ""; }
351 }
Int_t i
Definition: run_full.C:25
ArgVector fDaughterListNames
TH1F* PndListProvider::GetHisto ( )
inline

Definition at line 103 of file PndListProvider.h.

References fHisto.

103 {return fHisto;}
std::string PndListProvider::GetName ( )
inline

Definition at line 79 of file PndListProvider.h.

References fName.

79 {return fName;}
std::string fName
int PndListProvider::GetNCandidates ( )

Definition at line 214 of file PndListProvider.cxx.

References fOwnList, and RhoCandList::GetLength().

215 {
216  //if (fOwnList && fAntiOwnList) return (fOwnList->GetLength()+fAntiOwnList->GetLength());
217  //else
218  return fOwnList.GetLength();
219 }
Int_t GetLength() const
Definition: RhoCandList.h:46
RhoCandList fOwnList
int PndListProvider::GetNDaughters ( )
inline

Definition at line 98 of file PndListProvider.h.

References fNDaughters.

98 {return fNDaughters;}
int PndListProvider::GetType ( )

Definition at line 363 of file PndListProvider.cxx.

References fParticlePDG.

364 {
365  if (fParticlePDG) { return fParticlePDG->PdgCode(); }
366  else { return 0; }
367 }
TParticlePDG * fParticlePDG
bool PndListProvider::IsGeneric ( )
inline

Definition at line 80 of file PndListProvider.h.

References fIsGeneric.

Referenced by Print().

80 {return fIsGeneric;}
void PndListProvider::Print ( )

Definition at line 166 of file PndListProvider.cxx.

References fDaughterPDG, fDaughterPointers, fIsGeneric, fName, fNDaughters, fParticlePDG, i, and IsGeneric().

167 {
168  cout <<"[";
169  if (fParticlePDG) { cout <<fParticlePDG->PdgCode(); }
170  cout <<"]("<<fName<<")";
171  if (!fIsGeneric) {
172  cout <<" -> ";
173  for (int i=0; i<fNDaughters; i++) {
174  cout <<"[";
175  if (fDaughterPDG[i]) { cout << fDaughterPointers[i]->GetType(); }//cout <<fDaughterPDG[i]->PdgCode();
176  cout <<"]("<<fDaughterPointers[i]->GetName()<<")";//fDaughterListNames[i]<<") ";
177  if (fDaughterPointers[i]->IsGeneric()) { cout<<"G "; }
178  else { cout <<" "; }
179  }
180  }
181  cout <<endl;
182 }
Int_t i
Definition: run_full.C:25
std::vector< PndListProvider * > fDaughterPointers
std::string fName
std::vector< TParticlePDG * > fDaughterPDG
TParticlePDG * fParticlePDG
void PndListProvider::Reset ( )

Definition at line 186 of file PndListProvider.cxx.

References RhoCandList::Cleanup(), fCandsUptodate, fHistoFilled, and fOwnList.

Referenced by SetCandList(), and ~PndListProvider().

187 {
188  fOwnList.Cleanup();
189  fCandsUptodate=false;
190  fHistoFilled=false;
191 }
void Cleanup()
Definition: RhoCandList.cxx:62
RhoCandList fOwnList
void PndListProvider::SetCandList ( RhoCandList cl)

Definition at line 195 of file PndListProvider.cxx.

References RhoCandList::Append(), fCandsUptodate, fMinusSel, fOwnList, fParticlePDG, fPlusSel, and Reset().

196 {
197  Reset();
198 
199  if (fParticlePDG->AntiParticle()) {
200  if (fParticlePDG->Charge()<0) {
202  } else {
204  }
205  } else {
206  fOwnList=cl;
207  }
208 
209  fCandsUptodate=true;
210 }
RhoMinusParticleSelector * fMinusSel
void Append(const RhoCandidate *c)
Definition: RhoCandList.h:52
RhoCandList fOwnList
RhoPlusParticleSelector * fPlusSel
TParticlePDG * fParticlePDG
void PndListProvider::SetGeneric ( bool  isgeneric = true)
inline

Definition at line 85 of file PndListProvider.h.

References fIsGeneric.

85 {fIsGeneric=isgeneric;}
void PndListProvider::SetHisto ( TH1F *  h)
inline

Definition at line 87 of file PndListProvider.h.

References fHisto, and h.

87 {fHisto=h;}
void PndListProvider::SetMassSelector ( double  mean,
double  width 
)

Definition at line 381 of file PndListProvider.cxx.

References fMassSel, and fName.

382 {
383  if (!fMassSel) {
384  fMassSel=new RhoMassParticleSelector((fName+"sel").c_str(),mean,width);
385  }
386 }
RhoMassParticleSelector * fMassSel
std::string fName
Double_t mean[nsteps]
Definition: dedx_bands.C:65
void PndListProvider::SetToDump ( bool  todump = true)
inline

Definition at line 86 of file PndListProvider.h.

References fToDump.

86 {fToDump=todump;}
void PndListProvider::SetType ( std::string  pdgType)

Definition at line 120 of file PndListProvider.cxx.

References fParticlePDG.

121 {
122  fParticlePDG=TDatabasePDG::Instance()->GetParticle(pdgType.c_str());
123 }
TParticlePDG * fParticlePDG
void PndListProvider::SetType ( int  pdgcode)

Definition at line 127 of file PndListProvider.cxx.

References fParticlePDG.

128 {
129  fParticlePDG=TDatabasePDG::Instance()->GetParticle(pdgcode);
130 }
TParticlePDG * fParticlePDG
bool PndListProvider::ToDump ( )
inline

Definition at line 81 of file PndListProvider.h.

References fToDump.

81 {return fToDump;}

Member Data Documentation

bool PndListProvider::fCandsUptodate
private

Definition at line 130 of file PndListProvider.h.

Referenced by GetCandList(), PndListProvider(), Reset(), and SetCandList().

ArgVector PndListProvider::fDaughterListNames
private

Definition at line 123 of file PndListProvider.h.

Referenced by AddDecayProduct(), GetDecayProdName(), and PndListProvider().

std::vector<TParticlePDG*> PndListProvider::fDaughterPDG
private

Definition at line 122 of file PndListProvider.h.

Referenced by AddDaughterType(), GetDaughterType(), PndListProvider(), and Print().

std::vector<PndListProvider*> PndListProvider::fDaughterPointers
private

Definition at line 124 of file PndListProvider.h.

Referenced by AddDaughterPointer(), GetCandList(), and Print().

TH1F* PndListProvider::fHisto
private

Definition at line 119 of file PndListProvider.h.

Referenced by FillHisto(), GetCandList(), GetHisto(), PndListProvider(), and SetHisto().

bool PndListProvider::fHistoFilled
private

Definition at line 131 of file PndListProvider.h.

Referenced by FillHisto(), PndListProvider(), and Reset().

bool PndListProvider::fIsGeneric
private

Definition at line 128 of file PndListProvider.h.

Referenced by IsGeneric(), PndListProvider(), Print(), and SetGeneric().

RhoMassParticleSelector* PndListProvider::fMassSel
private

Definition at line 117 of file PndListProvider.h.

Referenced by GetCandList(), PndListProvider(), SetMassSelector(), and ~PndListProvider().

RhoMinusParticleSelector* PndListProvider::fMinusSel
private

Definition at line 115 of file PndListProvider.h.

Referenced by PndListProvider(), SetCandList(), and ~PndListProvider().

std::string PndListProvider::fName
private

Definition at line 111 of file PndListProvider.h.

Referenced by GetName(), PndListProvider(), Print(), and SetMassSelector().

int PndListProvider::fNDaughters
private
RhoCandList PndListProvider::fOwnList
private
TParticlePDG* PndListProvider::fParticlePDG
private

Definition at line 121 of file PndListProvider.h.

Referenced by GetType(), PndListProvider(), Print(), SetCandList(), and SetType().

RhoPlusParticleSelector* PndListProvider::fPlusSel
private

Definition at line 116 of file PndListProvider.h.

Referenced by PndListProvider(), SetCandList(), and ~PndListProvider().

bool PndListProvider::fToDump
private

Definition at line 129 of file PndListProvider.h.

Referenced by PndListProvider(), SetToDump(), and ToDump().


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