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

#include <PndSimpleCombiner.h>

Public Member Functions

 PndSimpleCombiner (PndAnalysis *fAna, TString decay, TString params="", double Ecm=0)
 
 ~PndSimpleCombiner ()
 
void SetVerbose (int verb=1)
 
void Combine ()
 
void Print ()
 
void SetPid (TString crit="", TString algo="")
 
void SetPidElectron (TString crit="", TString algo="")
 
void SetPidMuon (TString crit="", TString algo="")
 
void SetPidPion (TString crit="", TString algo="")
 
void SetPidKaon (TString crit="", TString algo="")
 
void SetPidProton (TString crit="", TString algo="")
 
int GetNLists ()
 
bool GetList (RhoCandList &l, TString comp)
 
bool GetList (RhoCandList &l, int pdg)
 
bool GetListN (RhoCandList &l, int idx)
 

Private Member Functions

void InitDecayInfo (SCDecayInfo &info, int pdg, int idx)
 
void FillGenericLists ()
 
int CombineList (RhoCandList &l, int mpdg, std::vector< int > &idx)
 
int SplitString (TString s, TString delim, StringList &toks)
 
bool ParseDecay (TString decay)
 
bool ParseParams (TString params)
 
bool CCInvariant (std::vector< int > &vpdg)
 
bool CCInvariant (int pdg)
 
int AntiPdg (int pdg)
 

Private Attributes

PndAnalysisfAnalysis
 
TString fDecay
 
TString fGlobParams
 
int fNLists
 
RhoCandList fList [MAXLISTS]
 
int fVerbose
 
double fEmin
 
double fPmin
 
double fEcm
 
RhoEnergyParticleSelectorfESel
 
RhoMomentumParticleSelectorfPSel
 
std::map< int, int > fPdgIdxMap
 
std::map< int, int > fIdxPdgMap
 
std::map< int, TStringfIdxListNameMap
 
std::map< int, TStringfIdxPidCritMap
 
std::map< int, TStringfIdxPidAlgoMap
 
std::vector< SCDecayInfofDecayInfoArray
 

Detailed Description

Definition at line 59 of file PndSimpleCombiner.h.

Constructor & Destructor Documentation

PndSimpleCombiner::PndSimpleCombiner ( PndAnalysis fAna,
TString  decay,
TString  params = "",
double  Ecm = 0 
)

Definition at line 55 of file PndSimpleCombiner.cxx.

References fEcm, fIdxListNameMap, fIdxPdgMap, fPdgIdxMap, i, ParseDecay(), ParseParams(), SetPid(), and TString.

55  :
56  fAnalysis(fAna), fDecay(decay), fGlobParams(params), fNLists(11), fVerbose(0),
57  fEmin(0.), fPmin(0.), fEcm(Ecm), fESel(0), fPSel(0)
58 {
59  // initialize mapping pdg -> list index and list name
60  int pdgcodes[] = { -11, 11, -13, 13, 211, -211, 321, -321, 2212, -2212, 22};
61  TString pdgnames[] = {"ElectronPlus", "ElectronMinus", "MuonPlus", "MuonMinus", "PionPlus", "PionMinus", "KaonPlus", "KaonMinus", "ProtonPlus", "ProtonMinus", "Neutral"};
62 
63  fPdgIdxMap.clear();
64  fIdxPdgMap.clear();
65  fIdxListNameMap.clear();
66 
67  for (int i=0;i<11;++i)
68  {
69  fPdgIdxMap[pdgcodes[i]] = i; // maps pdg code -> list index
70  fIdxPdgMap[i] = pdgcodes[i]; // maps list index -> pdg code
71  fIdxListNameMap[i] = pdgnames[i]; // maps list index -> generic list name (ElectronPlus, PionMinus, ...; see above)
72  }
73  // set initial pid configuration
74  SetPid("All", "PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts;PidAlgoSciT;PidAlgoRich;PidAlgoFtof");
75 
76  if (fEcm==0) fEcm=-999;
77 
78  assert(ParseDecay(decay));
79  ParseParams(params);
80 }
RhoMomentumParticleSelector * fPSel
PndAnalysis * fAnalysis
Int_t i
Definition: run_full.C:25
std::map< int, int > fIdxPdgMap
std::map< int, int > fPdgIdxMap
bool ParseDecay(TString decay)
RhoEnergyParticleSelector * fESel
void SetPid(TString crit="", TString algo="")
std::map< int, TString > fIdxListNameMap
bool ParseParams(TString params)
PndSimpleCombiner::~PndSimpleCombiner ( )

Definition at line 84 of file PndSimpleCombiner.cxx.

References fDecayInfoArray, fESel, fPSel, and i.

85 {
86  if (fESel) delete fESel;
87  if (fPSel) delete fPSel;
88 
89  for (size_t i=0;i<fDecayInfoArray.size();++i) if (fDecayInfoArray[i].msel) delete fDecayInfoArray[i].msel;
90 }
RhoMomentumParticleSelector * fPSel
Int_t i
Definition: run_full.C:25
std::vector< SCDecayInfo > fDecayInfoArray
RhoEnergyParticleSelector * fESel

Member Function Documentation

int PndSimpleCombiner::AntiPdg ( int  pdg)
private

Definition at line 475 of file PndSimpleCombiner.cxx.

Referenced by CCInvariant(), Combine(), FillGenericLists(), and ParseDecay().

476 {
477  if (TDatabasePDG::Instance()->GetParticle(pdg))
478  {
479  if (!TDatabasePDG::Instance()->GetParticle(pdg)->AntiParticle()) return pdg;
480  else return TDatabasePDG::Instance()->GetParticle(pdg)->AntiParticle()->PdgCode();
481  }
482  return -999999;
483 }
bool PndSimpleCombiner::CCInvariant ( std::vector< int > &  vpdg)
private

Definition at line 457 of file PndSimpleCombiner.cxx.

References i.

Referenced by ParseDecay().

458 {
459  int sum=0, nminus=0, nplus=0;
460 
461  for (size_t i=0;i<vpdg.size();++i)
462  {
463  if (!CCInvariant(vpdg[i]))
464  {
465  sum += vpdg[i];
466  if (vpdg[i]>0) nplus++;
467  else nminus++;
468  }
469  }
470  return (sum==0 && nplus==nminus);
471 }
Int_t i
Definition: run_full.C:25
bool CCInvariant(std::vector< int > &vpdg)
bool PndSimpleCombiner::CCInvariant ( int  pdg)
inlineprivate

Definition at line 92 of file PndSimpleCombiner.h.

References AntiPdg().

92 {return (pdg == AntiPdg(pdg));} // is pdg = pdg of anti-particle
void PndSimpleCombiner::Combine ( )

Definition at line 524 of file PndSimpleCombiner.cxx.

References AntiPdg(), RhoCandList::Append(), CombineList(), SCDecayInfo::daucc, SCDecayInfo::didx, SCDecayInfo::dpdg, fDecayInfoArray, FillGenericLists(), fList, fPdgIdxMap, i, SCDecayInfo::midx, SCDecayInfo::mpdg, SCDecayInfo::msel, n, and RhoCandList::Select().

Referenced by PndSimpleCombinerTask::Exec().

525 {
526  // fill list of final state particles
528 
529  // loop through list definitions and produce composites
530  int n = fDecayInfoArray.size();
531  for (int i=0;i<n;++i)
532  {
533  SCDecayInfo info = fDecayInfoArray[i];
534  CombineList(fList[info.midx], info.mpdg, info.didx);
535 
536  // do we need to add the cc FS?
537  if (info.daucc)
538  {
539  // create index list with cc daughters
540  std::vector<int> aidx;
541  for (size_t j=0;j<info.didx.size();++j)
542  {
543  int apdg = AntiPdg(info.dpdg[j]);
544  aidx.push_back(fPdgIdxMap[apdg]);
545  }
546  // create temporary list with the cc combinatorics
547  RhoCandList l;
548  CombineList(l, info.mpdg, aidx);
549  // append it to the original list
550  fList[info.midx].Append(l);
551  }
552  // if there is a mass selector connected, apply it
553  if (info.msel) fList[info.midx].Select(info.msel);
554  }
555 }
void Append(const RhoCandidate *c)
Definition: RhoCandList.h:52
Int_t i
Definition: run_full.C:25
std::vector< SCDecayInfo > fDecayInfoArray
int n
RhoCandList fList[MAXLISTS]
std::map< int, int > fPdgIdxMap
std::vector< int > dpdg
int CombineList(RhoCandList &l, int mpdg, std::vector< int > &idx)
void Select(RhoParticleSelectorBase *pidmgr)
std::vector< int > didx
RhoMassParticleSelector * msel
int PndSimpleCombiner::CombineList ( RhoCandList l,
int  mpdg,
std::vector< int > &  idx 
)
private

Definition at line 559 of file PndSimpleCombiner.cxx.

References RhoCandList::Cleanup(), RhoCandList::Combine(), fList, and RhoCandList::GetLength().

Referenced by Combine().

560 {
561  l.Cleanup();
562 
563  int nd = idx.size();
564 
565  switch (nd)
566  {
567  case 2: l.Combine(fList[idx[0]], fList[idx[1]], mpdg); break;
568  case 3: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], mpdg); break;
569  case 4: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], mpdg); break;
570  case 5: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], mpdg); break;
571 
572  case 6: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], fList[idx[5]], mpdg); break;
573  case 7: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], fList[idx[5]], fList[idx[6]], mpdg); break;
574  case 8: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], fList[idx[5]], fList[idx[6]], fList[idx[7]], mpdg); break;
575  case 9: l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], fList[idx[5]], fList[idx[6]], fList[idx[7]], fList[idx[8]], mpdg); break;
576  case 10:l.Combine(fList[idx[0]], fList[idx[1]], fList[idx[2]], fList[idx[3]], fList[idx[4]], fList[idx[5]], fList[idx[6]], fList[idx[7]], fList[idx[8]], fList[idx[9]], mpdg); break;
577  }
578 
579  return l.GetLength();
580 }
void Cleanup()
Definition: RhoCandList.cxx:62
Int_t GetLength() const
Definition: RhoCandList.h:46
RhoCandList fList[MAXLISTS]
int idx[MAX]
Definition: autocutx.C:38
void Combine(RhoCandList &l1, RhoCandList &l2)
void PndSimpleCombiner::FillGenericLists ( )
private

Definition at line 409 of file PndSimpleCombiner.cxx.

References AntiPdg(), SCDecayInfo::didx, SCDecayInfo::dpdg, fAnalysis, fDecayInfoArray, fESel, fIdxListNameMap, fIdxPdgMap, fIdxPidAlgoMap, fIdxPidCritMap, PndAnalysis::FillList(), fList, fPdgIdxMap, fPSel, fVerbose, RhoCandList::GetLength(), i, n, SCDecayInfo::ndaug, and RhoCandList::Select().

Referenced by Combine().

410 {
411  // fill all generic lists, which are used by any of the composites
412  // loop through composites
413  bool filled[11] = {0}; // marker, whether list is already filled
414 
415  int n = fDecayInfoArray.size();
416  for (int i=0;i<n;++i)
417  {
418  SCDecayInfo info = fDecayInfoArray[i];
419 
420  for (int j=0;j<info.ndaug;++j)
421  {
422  int lidx = info.didx[j];
423 
424  if (lidx>10 || filled[lidx]) continue; // this list is not generic or already filled
425 
426  filled[lidx] = true; // mark list as filled
427  int laidx = fPdgIdxMap[AntiPdg(info.dpdg[j])]; // index of antiparticle list (only relevant for charged)
428 
429  if (lidx<10) // charged lists
430  {
432  if (fPSel) fList[lidx].Select(fPSel);
433 
434  // also fill anti pdg generic list for charged generic ...
436  if (fPSel) fList[laidx].Select(fPSel);
437  filled[laidx] = true; // ... and mark list as filled
438  }
439  else if (lidx==10) //neutrals
440  {
441  fAnalysis->FillList(fList[lidx],fIdxListNameMap[lidx]);
442  if (fESel) fList[lidx].Select(fESel);
443  }
444  if (lidx<11 && fVerbose)
445  {
446  cout <<"idx :"<<lidx<<" pdg:"<<fIdxPdgMap[lidx]<<" N:"<<fList[lidx].GetLength()<<endl;
447  if (lidx<10) cout <<"aidx:"<<laidx<<" pdg:"<<fIdxPdgMap[laidx]<<" N:"<<fList[laidx].GetLength()<<endl;
448  }
449  }
450  }
451 }
RhoMomentumParticleSelector * fPSel
PndAnalysis * fAnalysis
std::map< int, TString > fIdxPidAlgoMap
Int_t i
Definition: run_full.C:25
std::map< int, TString > fIdxPidCritMap
Int_t GetLength() const
Definition: RhoCandList.h:46
std::vector< SCDecayInfo > fDecayInfoArray
int n
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
std::map< int, int > fIdxPdgMap
RhoCandList fList[MAXLISTS]
std::map< int, int > fPdgIdxMap
std::vector< int > dpdg
void Select(RhoParticleSelectorBase *pidmgr)
RhoEnergyParticleSelector * fESel
std::vector< int > didx
std::map< int, TString > fIdxListNameMap
bool PndSimpleCombiner::GetList ( RhoCandList l,
TString  comp 
)

Definition at line 585 of file PndSimpleCombiner.cxx.

Referenced by PndSimpleCombinerTask::Exec().

586 {
587  if (!TDatabasePDG::Instance()->GetParticle(comp)) return false;
588 
589  return GetList(l, TDatabasePDG::Instance()->GetParticle(comp)->PdgCode());
590 }
bool GetList(RhoCandList &l, TString comp)
bool PndSimpleCombiner::GetList ( RhoCandList l,
int  pdg 
)

Definition at line 594 of file PndSimpleCombiner.cxx.

References RhoCandList::Cleanup(), fList, and fPdgIdxMap.

595 {
596  l.Cleanup();
597 
598  if (fPdgIdxMap.find(pdg) == fPdgIdxMap.end()) return false;
599 
600  l = fList[fPdgIdxMap[pdg]];
601  return true;
602 }
void Cleanup()
Definition: RhoCandList.cxx:62
RhoCandList fList[MAXLISTS]
std::map< int, int > fPdgIdxMap
bool PndSimpleCombiner::GetListN ( RhoCandList l,
int  idx 
)

Definition at line 606 of file PndSimpleCombiner.cxx.

References RhoCandList::Cleanup(), fList, and GetNLists().

607 {
608  l.Cleanup();
609 
610  if (idx>=0 && idx<GetNLists())
611  {
612  l = fList[idx+11];
613  return true;
614  }
615 
616  return false;
617 }
void Cleanup()
Definition: RhoCandList.cxx:62
RhoCandList fList[MAXLISTS]
int idx[MAX]
Definition: autocutx.C:38
int PndSimpleCombiner::GetNLists ( )
inline

Definition at line 79 of file PndSimpleCombiner.h.

References fNLists.

Referenced by GetListN().

79 {return fNLists-11;}
void PndSimpleCombiner::InitDecayInfo ( SCDecayInfo info,
int  pdg,
int  idx 
)
private

Definition at line 114 of file PndSimpleCombiner.cxx.

References SCDecayInfo::daucc, SCDecayInfo::didx, SCDecayInfo::dpdg, idx, SCDecayInfo::midx, SCDecayInfo::mpdg, SCDecayInfo::msel, SCDecayInfo::mwin, SCDecayInfo::mwinhi, and SCDecayInfo::mwinlo.

Referenced by ParseDecay().

115 {
116  info.mpdg = pdg; // mother(composite) pdg code
117  info.midx = idx; // add another list for the new composite
118  info.daucc = false; // in case only the daughters have a cc, but not the mother, we have to care (e.g. etac -> Ks K+ pi-)
119  info.mwin = 0; // default: no mass selection
120  info.mwinlo = 0; // default: no mass selection
121  info.mwinhi = 0; // default: no mass selection
122  info.msel = 0; // default: no mass selection
123  info.dpdg.clear(); // pdgs of daughters
124  info.didx.clear(); // list index of daughters
125 }
std::vector< int > dpdg
int idx[MAX]
Definition: autocutx.C:38
std::vector< int > didx
RhoMassParticleSelector * msel
bool PndSimpleCombiner::ParseDecay ( TString  decay)
private

Definition at line 129 of file PndSimpleCombiner.cxx.

References AntiPdg(), CCInvariant(), SCDecayInfo::didx, SCDecayInfo::dpdg, fDecayInfoArray, fIdxPdgMap, fNLists, fPdgIdxMap, i, InitDecayInfo(), SCDecayInfo::midx, SCDecayInfo::mpdg, SCDecayInfo::ndaug, SplitString(), and TString.

Referenced by PndSimpleCombiner().

130 {
131  StringList subdec;
132 
133  int ndec = SplitString(decay, ";", subdec);
134 
135  // loop over subdecays
136  for (int i=0;i<ndec;++i)
137  {
138  StringList dectoks;
139  // does decay string contain exactly one '->' in the middle of the string?
140  subdec[i].ReplaceAll("->",">");
141  SplitString(subdec[i],">",dectoks);
142 
143  if ( dectoks.size()!=2 || !TDatabasePDG::Instance()->GetParticle(dectoks[0]) ) {cout <<"[PndSimpleCombiner] **** ERROR : Invalid decay pattern '"<<subdec[i].Data()<<"'"<<endl; return false;}
144 
145  // split string again; first token is supposed to be the decaying resonance
146  TString curstr = dectoks[0]+" "+dectoks[1];
147  SplitString(curstr," ",dectoks);
148 
149  // too many daughters (>10)
150  if (dectoks.size()>11) {cout <<"[PndSimpleCombiner] **** ERROR : Exceeding max number of daughters (10): '"<<subdec[i].Data()<<"'"<<endl; return false;}
151 
152  SCDecayInfo info;
153  InitDecayInfo(info, TDatabasePDG::Instance()->GetParticle(dectoks[0])->PdgCode(), fNLists++);
154 
155  fPdgIdxMap[info.mpdg] = info.midx; // add new list to the pdg <-> list index maps
156  fIdxPdgMap[info.midx] = info.mpdg;
157 
158  bool cc=true;
159 
160  // loop over daughters
161  for (size_t j=1;j<dectoks.size();++j)
162  {
163  if (dectoks[j]=="nocc") { cc = false; continue; }
164 
165  if (!TDatabasePDG::Instance()->GetParticle(dectoks[j])) {cout <<"[PndSimpleCombiner] **** ERROR : Unknown particle '"<<dectoks[j].Data()<<"'"<<endl; return false;}
166 
167  // pdg code of daughter
168  int dpdg = TDatabasePDG::Instance()->GetParticle(dectoks[j])->PdgCode();
169  // add corresponding index if exists
170  if ( fPdgIdxMap.find(dpdg) == fPdgIdxMap.end() ) {cout <<"[PndSimpleCombiner] **** ERROR : Undefined list '"<<dectoks[j].Data()<<"'"<<endl; return false;}
171  if ( dpdg == info.mpdg ) {cout <<"[PndSimpleCombiner] **** ERROR : Invalid recursion in '"<<subdec[i].Data()<<"'"<<endl; return false;}
172 
173  info.dpdg.push_back(dpdg);
174  info.didx.push_back(fPdgIdxMap[dpdg]);
175  }
176  info.ndaug = info.didx.size();
177  fDecayInfoArray.push_back(info);
178 
179  // do we need the cc decay definition
180  if (cc)
181  {
182  // does the cc mode exist?
183  // either the mother or the FS has to have a cc not being identical
184  // example: 1) D0 -> K- pi+ has cc -> D0bar -> K+ pi- (different FS and mother)
185  // 2) D0 -> KS pi+ pi- has cc -> D0bar -> KS pi+ pi- (different mother)
186  // 3) etac -> KS K- pi+ has cc -> etac -> KS K+ pi- (different FS)
187  // 4) phi -> K+ K- has _no_ cc! (identical FS and mother)
188 
189  // does the mother have a cc? (case 1 + 2)
190  // then we add a new list
191  if (!CCInvariant(info.mpdg))
192  {
193  SCDecayInfo ainfo;
194  InitDecayInfo(ainfo, AntiPdg(info.mpdg), fNLists++);
195 
196  fPdgIdxMap[ainfo.mpdg] = ainfo.midx;
197  fIdxPdgMap[ainfo.midx] = ainfo.mpdg;
198 
199  for (size_t j=0;j<info.dpdg.size();++j)
200  {
201  int apdg = AntiPdg(info.dpdg[j]);
202  if (apdg==-999999 || fPdgIdxMap.find(apdg) == fPdgIdxMap.end() ) {cout <<"[PndSimpleCombiner] **** ERROR : No list for PDG code "<<apdg<<endl; return false;}
203  ainfo.dpdg.push_back(apdg);
204  ainfo.didx.push_back(fPdgIdxMap[apdg]);
205  }
206  ainfo.ndaug = info.didx.size();
207  fDecayInfoArray.push_back(ainfo);
208  }
209  // only the daughters have a cc (case 3)
210  // then we simply set the daucc switch, which will be taken into account during combinatorics
211  else if (!CCInvariant(info.dpdg))
212  {
213  for (size_t j=0;j<info.dpdg.size();++j)
214  {
215  int apdg = AntiPdg(info.dpdg[j]);
216  if (apdg==-999999 || fPdgIdxMap.find(apdg) == fPdgIdxMap.end() ) {cout <<"[PndSimpleCombiner] **** ERROR : No list for PDG code "<<apdg<<endl; return false;}
217  }
218 
219  fDecayInfoArray[fDecayInfoArray.size()-1].daucc=true;
220  }
221  }
222  }
223 
224  return true;
225 }
Int_t i
Definition: run_full.C:25
std::vector< SCDecayInfo > fDecayInfoArray
std::map< int, int > fIdxPdgMap
std::map< int, int > fPdgIdxMap
int SplitString(TString s, TString delim, StringList &toks)
std::vector< int > dpdg
std::vector< TString > StringList
void InitDecayInfo(SCDecayInfo &info, int pdg, int idx)
bool CCInvariant(std::vector< int > &vpdg)
std::vector< int > didx
bool PndSimpleCombiner::ParseParams ( TString  params)
private

Definition at line 229 of file PndSimpleCombiner.cxx.

References fDecayInfoArray, fEcm, fEmin, fESel, fIdxListNameMap, fPdgIdxMap, fPmin, fPSel, i, mean, SCDecayInfo::mpdg, SCDecayInfo::msel, SCDecayInfo::mwin, SCDecayInfo::mwinhi, SCDecayInfo::mwinlo, SetPid(), SetPidElectron(), SetPidKaon(), SetPidMuon(), SetPidPion(), SetPidProton(), SplitString(), and TString.

Referenced by PndSimpleCombiner().

230 {
231  StringList parm;
232  SplitString(params,":",parm);
233 
234  for (size_t i=0;i<parm.size();++i)
235  {
236  if (parm[i]=="ebrem")
237  {
238  fIdxListNameMap[0]="ElectronPlusBrem";
239  fIdxListNameMap[1]="ElectronMinusBrem";
240  continue;
241  }
242 
243  StringList pair;
244  SplitString(parm[i],"=",pair);
245 
246  if (pair.size()!=2) {cout <<"[PndSimpleCombiner] **** WARNING : Invalid parameter setting '"<<parm[i]<<"' ignored"<<endl; continue;}
247 
248  // ****
249  // global mass window setting
250  // can be either : mwin(part)=0.5 -> cut +- 0.25 around nominal mass of part
251  // or : mwin(part)=0.25|0.75 -> cut 0.25 < mpart < 0.75
252  // ****
253  if (pair[0]=="mwin")
254  {
255  double window = pair[1].Atof();
256 
257  for (size_t j=0;j<fDecayInfoArray.size();++j)
258  {
259  SCDecayInfo &info = fDecayInfoArray[j];
260  if (info.msel) delete info.msel;
261 
262  // fetch center as nominal pdg mass; if pbarpSystemX (8888x) set center to Ecm
263  double mean = TDatabasePDG::Instance()->GetParticle(info.mpdg)->Mass();
264 
265  if (info.mpdg/10 == 8888) { mean = fEcm;}
266 
267  info.mwin = window;
268  info.mwinlo = mean - window/2;
269  info.mwinhi = mean + window/2;
270 
271  // check if Ecm was not set and particle is pbarpSystem
272  if (mean>-999)
273  info.msel = new RhoMassParticleSelector("msel",mean,window);
274  else
275  {
276  info.mwin = -1; // indicates that we skip the mass window cut for pbarpSystem if Ecm was not set
277  }
278  }
279  }
280  // ****
281  // mass window for one composite
282  // ****
283  else if (pair[0].BeginsWith("mwin"))
284  {
285  // extract particle name from strings like 'mwin(D0)'
286  pair[0] = pair[0](5,pair[0].Length()-6);
287 
288  if (!TDatabasePDG::Instance()->GetParticle(pair[0])) {cout <<"[PndSimpleCombiner] **** WARNING : Unknown particle type '"<<pair[0]<<"'"<<endl;continue;}
289 
290  int pdg = TDatabasePDG::Instance()->GetParticle(pair[0])->PdgCode();
291  if (fPdgIdxMap.find(pdg) == fPdgIdxMap.end()) {cout <<"[PndSimpleCombiner] **** WARNING : Unknown particle list '"<<pair[0]<<"'"<<endl;continue;}
292 
293  // set default mean and window
294  double mean = TDatabasePDG::Instance()->GetParticle(pdg)->Mass(), window = 0;
295  if (pdg/10 == 8888) { mean = fEcm; }
296 
297  // range setting like 1.2|1.6 ?
298  if (pair[1].Contains("|"))
299  {
300  double lower = TString((pair[1])(0,pair[1].Index("|"))).Atof();
301  double upper = TString((pair[1])(pair[1].Index("|")+1,1000)).Atof();
302 
303  mean = (upper+lower)/2.;
304  window = upper-lower;
305 
306  if (mean<0 || window <=0) {cout <<"[PndSimpleCombiner] **** WARNING : Invalid mass window defintion '"<<parm[i]<<"'"<<endl;continue;}
307  }
308  else window = pair[1].Atof();
309 
310  for (size_t j=0;j<fDecayInfoArray.size();++j)
311  {
312  SCDecayInfo &info = fDecayInfoArray[j];
313  // only set for this particle type
314  if (abs(info.mpdg) == abs(pdg))
315  {
316  if (info.msel) delete info.msel;
317  info.mwin = window;
318  info.mwinlo = mean - window/2;
319  info.mwinhi = mean + window/2;
320 
321  // check if Ecm was not set and particle is pbarpSystem
322  if (mean>-999)
323  info.msel = new RhoMassParticleSelector("msel",mean,window);
324  else
325  {
326  info.mwin = -1; // indicates that we skip the mass window cut for pbarpSystem if Ecm was not set
327  }
328  }
329  }
330  }
331 
332  // ****
333  // check for pid setting
334  // ****
335  if (pair[0] == "pid") SetPid(pair[1]);
336 
337  if (pair[0] == "pide") SetPidElectron(pair[1]);
338  if (pair[0] == "pidmu") SetPidMuon(pair[1]);
339  if (pair[0] == "pidpi") SetPidPion(pair[1]);
340  if (pair[0] == "pidk") SetPidKaon(pair[1]);
341  if (pair[0] == "pidp") SetPidProton(pair[1]);
342 
343  if (pair[0] == "algo") SetPid("",pair[1]);
344 
345  if (pair[0] == "algoe") SetPidElectron("",pair[1]);
346  if (pair[0] == "algomu") SetPidMuon("",pair[1]);
347  if (pair[0] == "algopi") SetPidPion("",pair[1]);
348  if (pair[0] == "algok") SetPidKaon("",pair[1]);
349  if (pair[0] == "algop") SetPidProton("",pair[1]);
350 
351  // ****
352  // E_min or p_min
353  // ****
354  if (pair[0] == "emin") {fEmin = pair[1].Atof(); fESel = new RhoEnergyParticleSelector("eSel",50.+fEmin,100.);}
355  if (pair[0] == "pmin") {fPmin = pair[1].Atof(); fPSel = new RhoMomentumParticleSelector("pSel",50.+fPmin,100.);}
356  }
357  return kTRUE;
358 }
RhoMomentumParticleSelector * fPSel
void SetPidPion(TString crit="", TString algo="")
void SetPidElectron(TString crit="", TString algo="")
Int_t i
Definition: run_full.C:25
void SetPidProton(TString crit="", TString algo="")
std::vector< SCDecayInfo > fDecayInfoArray
std::map< int, int > fPdgIdxMap
int SplitString(TString s, TString delim, StringList &toks)
std::vector< TString > StringList
RhoEnergyParticleSelector * fESel
Double_t mean[nsteps]
Definition: dedx_bands.C:65
void SetPidMuon(TString crit="", TString algo="")
RhoMassParticleSelector * msel
void SetPidKaon(TString crit="", TString algo="")
void SetPid(TString crit="", TString algo="")
std::map< int, TString > fIdxListNameMap
void PndSimpleCombiner::Print ( )

Definition at line 488 of file PndSimpleCombiner.cxx.

References SCDecayInfo::daucc, SCDecayInfo::didx, SCDecayInfo::dpdg, fDecayInfoArray, fEcm, fEmin, fIdxListNameMap, fIdxPidAlgoMap, fIdxPidCritMap, fPmin, i, SCDecayInfo::midx, SCDecayInfo::mpdg, SCDecayInfo::mwin, SCDecayInfo::mwinhi, SCDecayInfo::mwinlo, n, and printf().

Referenced by PndSimpleCombinerTask::Init().

489 {
490  cout <<"\n------------------------------------------"<<endl<<"[PndSimpleCombiner] **** Configuration"<<endl<<"------------------------------------------"<<endl;
491 
492  cout <<"E_cm = "<<fEcm<<endl<<endl;
493  cout <<"Neutrals E_min = "<<fEmin<<endl;
494  cout <<"Charged p_min = "<<fPmin<<endl<<endl;
495 
496  // fill all generic lists, which are used by any of the composites
497  // loop through composites
498  int n = fDecayInfoArray.size();
499  for (int i=0;i<n;++i)
500  {
501  SCDecayInfo info = fDecayInfoArray[i];
502  cout <<"Decay "<<i<<" : " << TDatabasePDG::Instance()->GetParticle(info.mpdg)->GetName()<<"("<<info.mpdg<<"/"<<info.midx<<") -> ";
503  for (size_t j=0;j<info.dpdg.size();++j) cout << TDatabasePDG::Instance()->GetParticle(info.dpdg[j])->GetName()<<"("<<info.dpdg[j]<<"/"<<info.didx[j]<<") ";
504 
505  cout <<" // daucc: "<<info.daucc;
506  cout <<" // mass window: ";
507  if (info.mwin>0) cout <<" ["<<info.mwinlo<<" ; "<<info.mwinhi<<"] = "<<info.mwin<<" GeV/cē"<<endl;
508  else if (info.mwin<0) cout <<" *** skipped due to missing Ecm setting *** "<<endl;
509  else cout <<"none"<<endl;
510  cout <<endl;
511  }
512  cout <<endl;
513  for (int i=0;i<10;++i)
514  {
515  printf("%-16s : %s, %s\n",fIdxListNameMap[i].Data(),fIdxPidCritMap[i].Data(),fIdxPidAlgoMap[i].Data());
516  }
517 
518  cout <<"\n------------------------------------------\n\n"<<endl;
519 }
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
std::map< int, TString > fIdxPidAlgoMap
Int_t i
Definition: run_full.C:25
std::map< int, TString > fIdxPidCritMap
std::vector< SCDecayInfo > fDecayInfoArray
int n
std::vector< int > dpdg
std::vector< int > didx
std::map< int, TString > fIdxListNameMap
void PndSimpleCombiner::SetPid ( TString  crit = "",
TString  algo = "" 
)

Definition at line 362 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, fIdxPidCritMap, SetPidElectron(), SetPidKaon(), SetPidMuon(), SetPidPion(), and SetPidProton().

Referenced by ParseParams(), and PndSimpleCombiner().

363 {
364  if (crit!="") fIdxPidCritMap.clear();
365  if (algo!="") fIdxPidAlgoMap.clear();
366 
367  SetPidElectron(crit, algo);
368  SetPidMuon(crit, algo);
369  SetPidPion(crit, algo);
370  SetPidKaon(crit, algo);
371  SetPidProton(crit, algo);
372 }
void SetPidPion(TString crit="", TString algo="")
void SetPidElectron(TString crit="", TString algo="")
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void SetPidProton(TString crit="", TString algo="")
void SetPidMuon(TString crit="", TString algo="")
void SetPidKaon(TString crit="", TString algo="")
void PndSimpleCombiner::SetPidElectron ( TString  crit = "",
TString  algo = "" 
)

Definition at line 376 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, and fIdxPidCritMap.

Referenced by ParseParams(), and SetPid().

377 {
378  if (crit!="") { fIdxPidCritMap[0] = crit; fIdxPidCritMap[1] = crit;}
379  if (algo!="") { fIdxPidAlgoMap[0] = algo; fIdxPidAlgoMap[1] = algo;}
380 }
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void PndSimpleCombiner::SetPidKaon ( TString  crit = "",
TString  algo = "" 
)

Definition at line 394 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, and fIdxPidCritMap.

Referenced by ParseParams(), and SetPid().

395 {
396  if (crit!="") {fIdxPidCritMap[7] = crit; fIdxPidCritMap[6] = crit;}
397  if (algo!="") {fIdxPidAlgoMap[7] = algo; fIdxPidAlgoMap[6] = algo;}
398 }
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void PndSimpleCombiner::SetPidMuon ( TString  crit = "",
TString  algo = "" 
)

Definition at line 382 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, and fIdxPidCritMap.

Referenced by ParseParams(), and SetPid().

383 {
384  if (crit!="") {fIdxPidCritMap[2] = crit; fIdxPidCritMap[3] = crit;}
385  if (algo!="") {fIdxPidAlgoMap[2] = algo; fIdxPidAlgoMap[3] = algo;}
386 }
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void PndSimpleCombiner::SetPidPion ( TString  crit = "",
TString  algo = "" 
)

Definition at line 388 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, and fIdxPidCritMap.

Referenced by ParseParams(), and SetPid().

389 {
390  if (crit!="") {fIdxPidCritMap[5] = crit; fIdxPidCritMap[4] = crit;}
391  if (algo!="") {fIdxPidAlgoMap[5] = algo; fIdxPidAlgoMap[4] = algo;}
392 }
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void PndSimpleCombiner::SetPidProton ( TString  crit = "",
TString  algo = "" 
)

Definition at line 400 of file PndSimpleCombiner.cxx.

References fIdxPidAlgoMap, and fIdxPidCritMap.

Referenced by ParseParams(), and SetPid().

401 {
402  if (crit!="") {fIdxPidCritMap[9] = crit; fIdxPidCritMap[8] = crit;}
403  if (algo!="") {fIdxPidAlgoMap[9] = algo; fIdxPidAlgoMap[8] = algo;}
404 }
std::map< int, TString > fIdxPidAlgoMap
std::map< int, TString > fIdxPidCritMap
void PndSimpleCombiner::SetVerbose ( int  verb = 1)
inline

Definition at line 69 of file PndSimpleCombiner.h.

References fVerbose.

Referenced by PndSimpleCombinerTask::Init().

69 {fVerbose = verb;}
int PndSimpleCombiner::SplitString ( TString  s,
TString  delim,
StringList toks 
)
private

Definition at line 94 of file PndSimpleCombiner.cxx.

References i, and TString.

Referenced by ParseDecay(), and ParseParams().

95 {
96  toks.clear();
97 
98  TObjArray *tok = s.Tokenize(delim);
99  int N = tok->GetEntries();
100 
101  for (int i=0;i<N;++i)
102  {
103  TString st = ((TObjString*)tok->At(i))->String();
104  st.ReplaceAll("\t","");
105  st = st.Strip(TString::kBoth);
106  if (st != "") toks.push_back(st);
107  }
108 
109  return toks.size();
110 }
Int_t i
Definition: run_full.C:25
TLorentzVector s
Definition: Pnd2DStar.C:50

Member Data Documentation

PndAnalysis* PndSimpleCombiner::fAnalysis
private

Definition at line 95 of file PndSimpleCombiner.h.

Referenced by FillGenericLists().

TString PndSimpleCombiner::fDecay
private

Definition at line 96 of file PndSimpleCombiner.h.

std::vector<SCDecayInfo> PndSimpleCombiner::fDecayInfoArray
private
double PndSimpleCombiner::fEcm
private

Definition at line 103 of file PndSimpleCombiner.h.

Referenced by ParseParams(), PndSimpleCombiner(), and Print().

double PndSimpleCombiner::fEmin
private

Definition at line 101 of file PndSimpleCombiner.h.

Referenced by ParseParams(), and Print().

RhoEnergyParticleSelector* PndSimpleCombiner::fESel
private

Definition at line 105 of file PndSimpleCombiner.h.

Referenced by FillGenericLists(), ParseParams(), and ~PndSimpleCombiner().

TString PndSimpleCombiner::fGlobParams
private

Definition at line 97 of file PndSimpleCombiner.h.

std::map<int,TString> PndSimpleCombiner::fIdxListNameMap
private

Definition at line 111 of file PndSimpleCombiner.h.

Referenced by FillGenericLists(), ParseParams(), PndSimpleCombiner(), and Print().

std::map<int,int> PndSimpleCombiner::fIdxPdgMap
private

Definition at line 109 of file PndSimpleCombiner.h.

Referenced by FillGenericLists(), ParseDecay(), and PndSimpleCombiner().

std::map<int,TString> PndSimpleCombiner::fIdxPidAlgoMap
private
std::map<int,TString> PndSimpleCombiner::fIdxPidCritMap
private
RhoCandList PndSimpleCombiner::fList[MAXLISTS]
private

Definition at line 99 of file PndSimpleCombiner.h.

Referenced by Combine(), CombineList(), FillGenericLists(), GetList(), and GetListN().

int PndSimpleCombiner::fNLists
private

Definition at line 98 of file PndSimpleCombiner.h.

Referenced by GetNLists(), and ParseDecay().

std::map<int,int> PndSimpleCombiner::fPdgIdxMap
private
double PndSimpleCombiner::fPmin
private

Definition at line 102 of file PndSimpleCombiner.h.

Referenced by ParseParams(), and Print().

RhoMomentumParticleSelector* PndSimpleCombiner::fPSel
private

Definition at line 106 of file PndSimpleCombiner.h.

Referenced by FillGenericLists(), ParseParams(), and ~PndSimpleCombiner().

int PndSimpleCombiner::fVerbose
private

Definition at line 100 of file PndSimpleCombiner.h.

Referenced by FillGenericLists(), and SetVerbose().


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