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

#include <PndSimpleCombinerTask.h>

Inheritance diagram for PndSimpleCombinerTask:

Public Member Functions

 PndSimpleCombinerTask (TString anadecay, TString anaparms, double p=0, int run=0, int mode=0)
 
 ~PndSimpleCombinerTask ()
 
void SetMultFactor (int fac)
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void SetPidAlgo (TString algo)
 
void SetVerbose (int v=1)
 

Private Member Functions

int SplitString (TString s, TString delim, std::vector< TString > &toks)
 
int CountChargedDaughters (RhoCandidate *c)
 
void InitParms ()
 
 ClassDef (PndSimpleCombinerTask, 1)
 

Private Attributes

int fVerbose
 
int fEvtCount
 
int fRun
 
int fMode
 
int fRunMult
 
TLorentzVector fIni
 
TString fAnaDecay
 
TString fAnaParms
 
int fNntp
 
TString fPidAlgo
 
bool fQaMC
 
bool fQaEventShape
 
int fQaRecoInfo
 
bool fQaEvShapeNtp
 
bool fFit4C
 
bool fBest4C
 
bool fFitVtx
 
double fFit4CChiCut
 
double fFitVtxChiCut
 
PndAnalysisfAnalysis
 
PndSimpleCombinerfSimpleCombiner
 
std::vector< int > vmpdg
 
std::vector< RhoTuple * > vntp
 
unsigned int fNodump
 
RhoTuplenmc
 
RhoTuplenevt
 
TClonesArray * fOnlineFilterInfo
 

Detailed Description

Definition at line 43 of file PndSimpleCombinerTask.h.

Constructor & Destructor Documentation

PndSimpleCombinerTask::PndSimpleCombinerTask ( TString  anadecay,
TString  anaparms,
double  p = 0,
int  run = 0,
int  mode = 0 
)

Definition at line 71 of file PndSimpleCombinerTask.cxx.

References fIni, mp, and sqrt().

71  :
72  FairTask("PndSimpleCombinerTask"), fVerbose(0), fEvtCount(0), fRun(run), fMode(mode), fRunMult(10000),
73  fAnaDecay(anadecay), fAnaParms(anaparms), fNntp(0),
74  fPidAlgo("PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts;PidAlgoSciT;PidAlgoRich;PidAlgoFtof"),
75  fQaMC(false), fQaEventShape(false), fQaRecoInfo(0), fQaEvShapeNtp(false), fFit4C(false),
76  fBest4C(false), fFitVtx(false), fFit4CChiCut(1e15), fFitVtxChiCut(1e8),
77  fAnalysis(0), fSimpleCombiner(0), fNodump(0), nmc(0), nevt(0)
78 {
79  fIni.SetXYZT(0,0,0,0);
80  double mp = 0.938272;
81  if (p>0.0001) fIni.SetXYZT(0,0,p, sqrt(p*p+mp*mp)+mp);
82 }
Int_t run
Definition: autocutx.C:47
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
static const double mp
Definition: mzparameters.h:11
Double_t p
Definition: anasim.C:58
Int_t mode
Definition: autocutx.C:47
PndSimpleCombiner * fSimpleCombiner
PndSimpleCombinerTask::~PndSimpleCombinerTask ( )

Definition at line 87 of file PndSimpleCombinerTask.cxx.

References fAnalysis, and fSimpleCombiner.

88 {
89  if (fSimpleCombiner) delete fSimpleCombiner;
90  delete fAnalysis;
91 }
PndSimpleCombiner * fSimpleCombiner

Member Function Documentation

PndSimpleCombinerTask::ClassDef ( PndSimpleCombinerTask  ,
 
)
private
int PndSimpleCombinerTask::CountChargedDaughters ( RhoCandidate c)
private

Definition at line 226 of file PndSimpleCombinerTask.cxx.

References RhoCandidate::Charge(), RhoCandidate::Daughter(), fabs(), i, and RhoCandidate::NDaughters().

Referenced by Exec().

227 {
228  int nd = 0;
229 
230  for (int i=0;i<c->NDaughters();++i)
231  {
232  if (fabs(c->Daughter(i)->Charge())>1e-6) nd++;
233  }
234 
235  return nd;
236 }
Int_t i
Definition: run_full.C:25
RhoCandidate * Daughter(Int_t n)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Int_t NDaughters() const
Double_t Charge() const
Definition: RhoCandidate.h:184
void PndSimpleCombinerTask::Exec ( Option_t *  opt)
virtual

Definition at line 261 of file PndSimpleCombinerTask.cxx.

References RhoKinFitter::Add4MomConstraint(), all, RhoTuple::Column(), PndSimpleCombiner::Combine(), CountChargedDaughters(), RhoCandidate::Daughter(), RhoTuple::DumpData(), fAnalysis, fBest4C, fEvtCount, fFit4C, fFit4CChiCut, fFitVtx, fFitVtxChiCut, PndAnalysis::FillList(), fIni, RhoKinFitter::Fit(), RhoFitterBase::Fit(), fMode, fNntp, fNodump, fOnlineFilterInfo, fPidAlgo, fQaEventShape, fQaEvShapeNtp, fQaRecoInfo, fRun, fRunMult, fSimpleCombiner, RhoFitterBase::GetChi2(), PndAnalysis::GetEventInTask(), RhoCandidate::GetFit(), RhoCandList::GetLength(), PndSimpleCombiner::GetList(), RhoCandidate::GetMcTruth(), PndOnlineFilterInfo::GetNTagTotal(), i, mclist, mp, RhoCandidate::NDaughters(), nevt, nmc, p, RhoCandidate::P4(), sqrt(), PndOnlineFilterInfo::Tagged(), vmpdg, and vntp.

262 {
263  // *** some variables
264  int i=0,j=0;
265 
266  // *** necessary to read the next event
268 
269  // *** print event counter
270  if (!(++fEvtCount%100)) cout << "[PndSimpleCombinerTask] evt "<<fEvtCount<<endl;
271 
272  // *******
273  // ******* PUT ANALYSIS CODE HERE
274  // *******
275 
276  PndOnlineFilterInfo *stInfo = 0;
277  if (fOnlineFilterInfo) stInfo = ( PndOnlineFilterInfo* ) fOnlineFilterInfo->At ( 0 );
278 
279  // determine pbarmom if not done so far
280  if (fIni.T()<1e-6)
281  {
283  fAnalysis->FillList(mclist, "McTruth");
284  if (mclist.GetLength()>0)
285  {
286  double p = mclist[0]->P3().Pz();
287  double mp = 0.938272;
288  fIni.SetXYZT(0,0,p, sqrt(p*p+mp*mp)+mp);
289  }
290  }
291 
292  // *** QA tool for simple dumping of analysis results in RhoRuple
293  // *** WIKI: https://panda-wiki.gsi.de/foswiki/bin/view/Computing/PandaRootAnalysisJuly13#PndRhoTupleQA
294  PndRhoTupleQA qa(fAnalysis,fIni.P());
295 
296  // *** RhoCandLists for the analysis
297  RhoCandList l1, l2, all;
298 
299  // *** store MC info in ntuple
300  if (nmc)
301  {
302  nmc->Column("ev", (Int_t) fEvtCount);
303  qa.qaMcList(nmc);
304  nmc->DumpData();
305  }
306 
307  // *** Setup event shape object
308  PndEventShape *evsh = 0;
309 
311  {
312  fAnalysis->FillList(all, "All", fPidAlgo);
313  evsh = new PndEventShape(all, fIni, 0.05, 0.1);
314  }
315 
316  // *** Store event shape info in ntuple
317  if (fQaEvShapeNtp)
318  {
319  nevt->Column("ev", (Int_t) fEvtCount);
320  nevt->Column("run", (Int_t) fRun);
321  nevt->Column("uid", (Int_t) fRun*fRunMult+fEvtCount);
322  nevt->Column("mode", (Int_t) fMode);
323 
324  qa.qaP4("beam", fIni, nevt);
325  qa.qaEventShape("es", evsh, nevt);
326 
327  nevt->DumpData();
328  }
329 
330  // *****
331  // *** combinatorics
332  // *****
333 
334  if (fSimpleCombiner)
335  {
337 
338  // ntuple dump
339  for (i=0;i<fNntp;++i)
340  {
341  if (fNodump & (1<<i)) continue; // if dump of ntuple 'ntpi' is skipped continue
342 
343  int pdg = vmpdg[i];
344  int apdg = 0;
345  if (TDatabasePDG::Instance()->GetParticle(pdg)->AntiParticle()) apdg = TDatabasePDG::Instance()->GetParticle(pdg)->AntiParticle()->PdgCode();
346 
347  // check whether there is an own ntuple connected to the anti-particle pdg; if yes, reset apdg
348  for (j=0; j<fNntp; ++j) if (vmpdg[j]==apdg) {apdg=0; j=fNntp+1;}
349 
350  // merge list from particles and anti-particles
351  fSimpleCombiner->GetList(l1, pdg);
352  if (apdg!=0 && fSimpleCombiner->GetList(l2, apdg)) l1.Append(l2);
353 
354  //RhoMassParticleSelector msel("msel",TDatabasePDG::Instance()->GetParticle(pdg)->Mass(),0.2);
355  //l1.Select(&msel);
356 
357  // number of charged daughters for vtx fit
358  int ncdau = -1;
359 
360  // determine the best chi2 from 4C fit in case we only want to store the best candidate
361  double best4cChi2 = 1e10;
362 
363  for (j=0;j<l1.GetLength();++j)
364  {
365  if (ncdau<0) ncdau = CountChargedDaughters(l1[j]);
366 
367  Float_t mmiss = (fIni-(l1[j]->P4())).M();
368  Float_t msum = l1[j]->M() + mmiss;
369 
370  // in case we do 4C fit, we need to do determine first, whether the candidate is the best one
371  // flag whether the fit is accepted
372  bool fitaccept = true;
373 
374  // for the last list we perform a 4C fit
375  if (fFit4C && i==fNntp-1)
376  {
377  RhoKinFitter fit4c(l1[j]);
378  fit4c.Add4MomConstraint(fIni);
379  fit4c.Fit();
380 
381  double chi2_4c = fit4c.GetChi2();
382 
383  if (chi2_4c>=fFit4CChiCut) fitaccept = false;
384  if (chi2_4c>best4cChi2 || !fitaccept) continue;
385 
386  best4cChi2 = chi2_4c;
387 
388  RhoCandidate *cfit = l1[j]->GetFit();
389 
390  vntp[i]->Column("chi24c", (Float_t) chi2_4c);
391  qa.qaP4("f4cx", cfit->P4(), vntp[i]);
392 
393  for (int k=0;k<cfit->NDaughters();++k)
394  {
395  RhoCandidate *d0fit = cfit->Daughter(k);
396  qa.qaP4(TString::Format("f4cxd%d",k),d0fit->P4(),vntp[i]);
397 
398  for (int k2=0;k2<d0fit->NDaughters();++k2)
399  {
400  RhoCandidate *ddfit = d0fit->Daughter(k2);
401  qa.qaP4(TString::Format("f4cxd%dd%d",k,k2),ddfit->P4(),vntp[i]);
402  }
403  }
404  }
405 
406  vntp[i]->Column("ev", (Int_t) fEvtCount);
407  vntp[i]->Column("cand", (Int_t) j);
408  vntp[i]->Column("ncand", (Int_t) l1.GetLength());
409  vntp[i]->Column("run", (Int_t) fRun);
410  vntp[i]->Column("uid", (Int_t) fRun*fRunMult+fEvtCount);
411  vntp[i]->Column("mode", (Int_t) fMode);
412 
413  vntp[i]->Column("mmiss", (Float_t) mmiss);
414  vntp[i]->Column("msum", (Float_t) msum);
415 
416  qa.qaP4("beam", fIni, vntp[i]);
417 
418  // store information about composite candidate tree recursively (see analysis/AnalysisTools/PndRhoTupleQA)
419  qa.qaComp("x", l1[j], vntp[i]);
420  // also optional for reco info
421  if (fQaRecoInfo==1) qa.qaRecoShortTree("x", l1[j], vntp[i]);
422  else if (fQaRecoInfo==2) qa.qaRecoFullTree("x", l1[j], vntp[i]);
423 
424  // store info about event shapes
425  if (fQaEventShape) qa.qaEventShapeShort("es",evsh, vntp[i]);
426 
427  // *** store info from trigger
428  if (stInfo)
429  {
430  vntp[i]->Column("trig", (Int_t) stInfo->Tagged() ); // event triggered
431  vntp[i]->Column("ntrig", (Int_t) stInfo->GetNTagTotal()); // total number of triggered candidates from all active lines
432  }
433 
434  // store the 4-vector of the truth matched candidate (or a dummy, if not matched to keep ntuple consistent)
435  RhoCandidate *truth = l1[j]->GetMcTruth();
436  TLorentzVector lv;
437  if (truth) lv = truth->P4();
438  qa.qaP4("trx", lv, vntp[i]);
439 
440 
441  // shall we do a vertex fit?
442  if (fFitVtx && ncdau>1)
443  {
444  RhoKinVtxFitter vtxfitter(l1[j]); // *** instantiate the vertex fitter; input is the object to be fitted
445  vtxfitter.Fit(); // *** perform fit
446 
447  RhoCandidate *cfit = l1[j]->GetFit(); // *** get the fitted candidate
448 
449  qa.qaVtx("fvxx",cfit,vntp[i]);
450  qa.qaP4("fvxx", cfit->P4(), vntp[i]);
451  double chi2_vtx = vtxfitter.GetChi2(); // *** and the chi^2 of the fit
452  vntp[i]->Column("chi2vx", (Float_t) chi2_vtx);
453 
454  if (chi2_vtx>=fFitVtxChiCut) fitaccept = false;
455  }
456 
457  if (fitaccept && ((i<fNntp-1) || !fBest4C )) vntp[i]->DumpData();
458  }
459 
460  if (fBest4C && i==fNntp-1 && best4cChi2<1e10) { vntp[i]->DumpData();}
461  }
462  }
463 
464  if (evsh) delete evsh;
465 }
CandList mclist
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
RhoCandidate * Daughter(Int_t n)
static const double mp
Definition: mzparameters.h:11
Double_t p
Definition: anasim.C:58
int CountChargedDaughters(RhoCandidate *c)
void GetEventInTask()
TLorentzVector P4() const
Definition: RhoCandidate.h:195
void Column(const char *label, Bool_t value, Bool_t defval=0, const char *block=0)
Definition: RhoTuple.cxx:56
RhoCandidate * GetMcTruth() const
Definition: RhoCandidate.h:437
void DumpData()
Definition: RhoTuple.cxx:391
PndSimpleCombiner * fSimpleCombiner
Int_t NDaughters() const
RhoCandidate * GetFit() const
Definition: RhoCandidate.h:293
bool GetList(RhoCandList &l, TString comp)
std::vector< RhoTuple * > vntp
void PndSimpleCombinerTask::Finish ( )
virtual

Definition at line 468 of file PndSimpleCombinerTask.cxx.

References fNntp, fNodump, RhoTuple::GetInternalTree(), i, nevt, nmc, and vntp.

469 {
470  if (nmc) nmc->GetInternalTree()->Write();
471  if (nevt) nevt->GetInternalTree()->Write();
472 
473  for (int i=0;i<fNntp;++i) if (!(fNodump & (1<<i))) vntp[i]->GetInternalTree()->Write();
474 }
Int_t i
Definition: run_full.C:25
TTree * GetInternalTree()
Definition: RhoTuple.h:207
std::vector< RhoTuple * > vntp
InitStatus PndSimpleCombinerTask::Init ( )
virtual

Definition at line 146 of file PndSimpleCombinerTask.cxx.

References fAnaDecay, fAnalysis, fAnaParms, fFit4C, fFit4CChiCut, fFitVtx, fFitVtxChiCut, fIni, fNntp, fNodump, fOnlineFilterInfo, fQaEvShapeNtp, fQaMC, fSimpleCombiner, fVerbose, RhoTuple::GetInternalTree(), i, InitParms(), n, nevt, nmc, PndSimpleCombiner::Print(), PndSimpleCombiner::SetVerbose(), SplitString(), TString, vmpdg, and vntp.

147 {
148  InitParms();
149 
150  // *** initialize PndAnalysis object and SimpleCombiner
151  fAnalysis = new PndAnalysis();
152 
153  if (fAnaDecay!="")
154  {
156 
159  }
160 
161  // *******
162  // ******* PREPARE/CREATE THE STUFF YOU NEED
163  // *******
164 
165  // ***
166  // *** Prepare RhoTuple output
167  // ***
168 
169  // *** Save current gDirectory
170  TDirectory *dir = gDirectory;
171  FairRootManager::Instance()->GetOutFile()->cd();
172 
173  std::vector<TString> toks;
174  fNntp = SplitString(fAnaDecay,";",toks);
175 
176  // *** create some ntuples
177  for (int i=0;i<fNntp;++i)
178  {
179  RhoTuple *n = 0;
180  if (!(fNodump & (1<<i))) // do we write this ntuple?
181  {
182  n = new RhoTuple(TString::Format("ntp%d",i), toks[i]);
183  n->GetInternalTree()->SetDirectory(gDirectory);
184  }
185  vntp.push_back(n);
186 
187  TString pname = toks[i](0,toks[i].Index("->"));
188  pname.ReplaceAll(" ","");
189  if (TDatabasePDG::Instance()->GetParticle(pname)) {vmpdg.push_back(TDatabasePDG::Instance()->GetParticle(pname)->PdgCode());}
190  }
191 
192  // *** create MC ntuple
193  if (fQaMC) nmc = new RhoTuple("nmc", "mctruth info");
194  if (nmc) nmc->GetInternalTree()->SetDirectory(gDirectory);
195 
196  // *** create EventShape ntuple
197  if (fQaEvShapeNtp) nevt = new RhoTuple("nevt", "event shape info");
198  if (nevt) nevt->GetInternalTree()->SetDirectory(gDirectory);
199 
200  // *** restore original gDirectory
201  dir->cd();
202 
203  // *** Connect to the Online Filter Info
204  fOnlineFilterInfo = ( TClonesArray* ) FairRootManager::Instance()->GetObject ( "OnlineFilterInfo" );
205 
206 
207  // ****** Print out some info from PndSimpleCombinerTask
208  cout <<"\n------------------------------------------"<<endl<<"[PndSimpleCombinerTask] **** Configuration"<<endl<<"------------------------------------------"<<endl;
209  cout <<"Fitting : ";
210  if (fFit4C) cout <<"4-C ( chi^2 < "<<fFit4CChiCut<<")";
211  if (fFitVtx) cout <<" Vertex ( chi^2 < "<<fFitVtxChiCut<<")";
212  cout <<endl;
213 
214  cout <<"Ntuple output : ";
215  if (fQaMC) cout <<"nmc ";
216  if (fQaEvShapeNtp) cout <<"nevt ";
217  for (int i=0;i<fNntp;++i) if (!(fNodump & (1<<i))) cout <<"ntp"<<i<<"("<<TDatabasePDG::Instance()->GetParticle(vmpdg[i])->GetName()<<") ";
218  cout <<"\n------------------------------------------\n\n"<<endl;
219 
220 
221  return kSUCCESS;
222 }
int SplitString(TString s, TString delim, std::vector< TString > &toks)
Int_t i
Definition: run_full.C:25
int n
PndSimpleCombiner * fSimpleCombiner
TTree * GetInternalTree()
Definition: RhoTuple.h:207
void SetVerbose(int verb=1)
std::vector< RhoTuple * > vntp
void PndSimpleCombinerTask::InitParms ( )
private

Definition at line 93 of file PndSimpleCombinerTask.cxx.

References cut, fAnaParms, fBest4C, fFit4C, fFit4CChiCut, fFitVtx, fFitVtxChiCut, fNodump, fQaEventShape, fQaEvShapeNtp, fQaMC, fQaRecoInfo, i, SplitString(), and TString.

Referenced by Init().

94 {
95  fAnaParms.ReplaceAll(" ","");
96 
97  StringList pars;
98  SplitString(fAnaParms,":",pars);
99 
100  // loop over all parameters
101  for (unsigned int i=0;i<pars.size();++i)
102  {
103  // is this a parameter handled by task or by combiner? If yes, it will be deleted later from the list
104  bool taskparm = false;
105 
106  // simple parameter flag
107  if (pars[i]=="qamc") { fQaMC = true; taskparm = true; } // write mc information
108  if (pars[i]=="qaevtshape") { fQaEventShape = true; taskparm = true; } // write event shape info
109  if (pars[i]=="qaevs") { fQaEventShape = true; taskparm = true; } // write event shape info
110  if (pars[i]=="qarec") { fQaRecoInfo = 1 ; taskparm = true; } // write reco info summary (EMC, STT, ...)
111  if (pars[i]=="qarecfull") { fQaRecoInfo = 2 ; taskparm = true; } // write full reco info
112 
113  if (pars[i]=="nevt") { fQaEvShapeNtp = true; taskparm = true; } // write event shape info to extra ntuple
114  if (pars[i].Contains("fit4c")) { fFit4C = true; taskparm = true; } // perform 4c fit for last particle (usually pbarpSystemX)
115  if (pars[i].Contains("fit4cbest")) { fBest4C = true; taskparm = true; } // only store best 4C fitted candidate
116  if (pars[i].Contains("fitvtx")) { fFitVtx = true; taskparm = true; } // perform vertex fit if possible (not cascaded)
117 
118  if (pars[i].BeginsWith("!ntp")) // !ntpX avoids dump of ntuple of X-th resonance to save disc space (for e.g. subresonances of decay tree)
119  {
120  int ntpnum = ((TString)(pars[i](4,100))).Atoi();
121  if (ntpnum>=0 && ntpnum<32) fNodump |= (1<<ntpnum); // this is a bit marker; if i-th bit set, dump of ntpi is skipped
122  taskparm = true;
123  }
124 
125  // parameter pair for chi2 cut from fit?
126  if (pars[i].Contains("<"))
127  {
128  // extract the value
129  double cut = TString((pars[i])(pars[i].Index("<")+1,1000)).Atof();
130  if (cut<=0) { cout <<"[PndSimpleCombinerTask] **** ERROR : Invalid parameter setting '"<<pars[i]<<"'"<<endl;continue; }
131 
132  // set the cut for the corresponding fitter
133  if (pars[i].Contains("fit4c")) fFit4CChiCut=cut;
134  if (pars[i].Contains("fitvtx")) fFitVtxChiCut=cut;
135  }
136 
137  // delete parameter from list for PndSimpleCombiner
138  if (taskparm) fAnaParms.ReplaceAll(pars[i],"");
139  }
140 }
int SplitString(TString s, TString delim, std::vector< TString > &toks)
Int_t i
Definition: run_full.C:25
double cut[MAX]
Definition: autocutx.C:36
std::vector< TString > StringList
void PndSimpleCombinerTask::SetMultFactor ( int  fac)
inline

Definition at line 54 of file PndSimpleCombinerTask.h.

References fRunMult.

54 {fRunMult=fac;} // set run multiplicator
void PndSimpleCombinerTask::SetPidAlgo ( TString  algo)
inline

Definition at line 64 of file PndSimpleCombinerTask.h.

References fPidAlgo.

Referenced by prod_ana(), quickana(), and quickfsimana().

64 { fPidAlgo = algo;}
void PndSimpleCombinerTask::SetVerbose ( int  v = 1)
inline

Definition at line 65 of file PndSimpleCombinerTask.h.

References fVerbose, and v.

65 {fVerbose = v;} // verbosity level, also passed to PndSimpleCombiner
__m128 v
Definition: P4_F32vec4.h:4
int PndSimpleCombinerTask::SplitString ( TString  s,
TString  delim,
std::vector< TString > &  toks 
)
private

Definition at line 240 of file PndSimpleCombinerTask.cxx.

References i, and TString.

Referenced by Init(), and InitParms().

241 {
242  toks.clear();
243 
244  TObjArray *tok = s.Tokenize(delim);
245  int N = tok->GetEntries();
246 
247  for (int i=0;i<N;++i)
248  {
249  TString st = ((TObjString*)tok->At(i))->String();
250  st.ReplaceAll("\t","");
251  st = st.Strip(TString::kBoth);
252  if (st != "") toks.push_back(st);
253  }
254 
255  return toks.size();
256 }
Int_t i
Definition: run_full.C:25
TLorentzVector s
Definition: Pnd2DStar.C:50

Member Data Documentation

TString PndSimpleCombinerTask::fAnaDecay
private

Definition at line 82 of file PndSimpleCombinerTask.h.

Referenced by Init().

PndAnalysis* PndSimpleCombinerTask::fAnalysis
private

Definition at line 97 of file PndSimpleCombinerTask.h.

Referenced by Exec(), Init(), and ~PndSimpleCombinerTask().

TString PndSimpleCombinerTask::fAnaParms
private

Definition at line 83 of file PndSimpleCombinerTask.h.

Referenced by Init(), and InitParms().

bool PndSimpleCombinerTask::fBest4C
private

Definition at line 91 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and InitParms().

int PndSimpleCombinerTask::fEvtCount
private

Definition at line 76 of file PndSimpleCombinerTask.h.

Referenced by Exec().

bool PndSimpleCombinerTask::fFit4C
private

Definition at line 90 of file PndSimpleCombinerTask.h.

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

double PndSimpleCombinerTask::fFit4CChiCut
private

Definition at line 93 of file PndSimpleCombinerTask.h.

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

bool PndSimpleCombinerTask::fFitVtx
private

Definition at line 92 of file PndSimpleCombinerTask.h.

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

double PndSimpleCombinerTask::fFitVtxChiCut
private

Definition at line 94 of file PndSimpleCombinerTask.h.

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

TLorentzVector PndSimpleCombinerTask::fIni
private

Definition at line 81 of file PndSimpleCombinerTask.h.

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

int PndSimpleCombinerTask::fMode
private

Definition at line 78 of file PndSimpleCombinerTask.h.

Referenced by Exec().

int PndSimpleCombinerTask::fNntp
private

Definition at line 84 of file PndSimpleCombinerTask.h.

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

unsigned int PndSimpleCombinerTask::fNodump
private

Definition at line 102 of file PndSimpleCombinerTask.h.

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

TClonesArray* PndSimpleCombinerTask::fOnlineFilterInfo
private

Definition at line 106 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and Init().

TString PndSimpleCombinerTask::fPidAlgo
private

Definition at line 85 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and SetPidAlgo().

bool PndSimpleCombinerTask::fQaEventShape
private

Definition at line 87 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and InitParms().

bool PndSimpleCombinerTask::fQaEvShapeNtp
private

Definition at line 89 of file PndSimpleCombinerTask.h.

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

bool PndSimpleCombinerTask::fQaMC
private

Definition at line 86 of file PndSimpleCombinerTask.h.

Referenced by Init(), and InitParms().

int PndSimpleCombinerTask::fQaRecoInfo
private

Definition at line 88 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and InitParms().

int PndSimpleCombinerTask::fRun
private

Definition at line 77 of file PndSimpleCombinerTask.h.

Referenced by Exec().

int PndSimpleCombinerTask::fRunMult
private

Definition at line 79 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and SetMultFactor().

PndSimpleCombiner* PndSimpleCombinerTask::fSimpleCombiner
private

Definition at line 98 of file PndSimpleCombinerTask.h.

Referenced by Exec(), Init(), and ~PndSimpleCombinerTask().

int PndSimpleCombinerTask::fVerbose
private

Definition at line 75 of file PndSimpleCombinerTask.h.

Referenced by Init(), and SetVerbose().

RhoTuple* PndSimpleCombinerTask::nevt
private

Definition at line 104 of file PndSimpleCombinerTask.h.

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

RhoTuple* PndSimpleCombinerTask::nmc
private

Definition at line 103 of file PndSimpleCombinerTask.h.

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

std::vector<int> PndSimpleCombinerTask::vmpdg
private

Definition at line 100 of file PndSimpleCombinerTask.h.

Referenced by Exec(), and Init().

std::vector<RhoTuple*> PndSimpleCombinerTask::vntp
private

Definition at line 101 of file PndSimpleCombinerTask.h.

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


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