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

#include <PndFastSim.h>

Inheritance diagram for PndFastSim:

Public Member Functions

 PndFastSim (bool persist=true)
 
 ~PndFastSim ()
 
void Register ()
 
virtual InitStatus Init ()
 
virtual void Finish ()
 
virtual void Exec (Option_t *opt)
 
bool AddDetector (std::string name, std::string params="")
 
bool AddDetector (PndFsmAbsDet *det)
 
void SetVerbosity (int vb)
 
bool EnableSplitoffs (std::string fname="splitpars.dat")
 
bool MergeNeutralClusters (bool merge=true, double par=0.389)
 
void EnableElectronBremsstrahlung (bool brems=true)
 
void EnablePropagation (bool propagate=true, bool tostartvtx=true, bool usecovmatrix=true, double tolerance=0.0)
 
void SetUseFlatCov (bool v=true)
 
void SetMultFilter (TString type, int min, int max=1000)
 
void SetInvMassFilter (TString filter, double min, double max, int mult=1)
 
void SetSeed (unsigned int seed=65539)
 

Private Member Functions

PndFsmResponsesumResponse (FsmResponseList respList)
 
int acceptFilters (RhoCandList &l)
 
int chCon (int i)
 
void copyAndSetMass (RhoCandList &l, RhoCandList &nl, int pdg)
 
bool cutAndSmear (PndFsmTrack *t, PndFsmResponse *r)
 
bool cutAndSmear (PndFsmTrack *t)
 
void smearEnergy (PndFsmTrack *t, double dE)
 
void smearMomentum (PndFsmTrack *t, double dp)
 
void smearTheta (PndFsmTrack *t, double dtheta)
 
void smearPhi (PndFsmTrack *t, double dphi)
 
void smearVertex (PndFsmTrack *t, TVector3 dV)
 
void smearM (PndFsmTrack *t, double dm)
 
void smearM2 (PndFsmTrack *t, double m2)
 
void smearMvddEdx (PndFsmTrack *t, double dedx)
 
void smearTpcdEdx (PndFsmTrack *t, double dedx)
 
void smearSttdEdx (PndFsmTrack *t, double dedx)
 
void SetFlatCovMatrix (PndFsmTrack *t, double dp=0., double dtheta=0., double dphi=0., double dE=0., double dx=0., double dy=0., double dz=0.)
 
void UpdateGammaHit (PndFsmTrack *t)
 
virtual void SetParContainers ()
 
bool smearTrack (PndFsmTrack *t, int idx=-1)
 
 ClassDef (PndFastSim, 1)
 

Private Attributes

TClonesArray * fMcCandidates
 
TClonesArray * fPidChargedCand
 
TClonesArray * fPidNeutralCand
 
TClonesArray * fMicroCandidates
 
TClonesArray * fPidChargedProb
 
TClonesArray * fPidNeutralProb
 PndPidProbability TCA for charged particles. More...
 
std::map< TString, TClonesArray * > fPidArrayList
 PndPidProbability TCA for neutral particles. More...
 
TClonesArray * fEventInfo
 PndPidProbability TCA's for individual detectors. More...
 
TRandom3 * fRand
 
int fVb
 
int evtcnt
 
TF1 * fspo [5][4]
 
TF1 * fBremsEnergy
 
bool fGenSplitOffs
 
bool fMergeNeutralClusters
 
bool fElectronBrems
 
double fMergeProbPar
 
bool fPropagate
 
bool fToStartVtx
 
bool fUseCovMatrix
 
bool fUseFlatCovMatrix
 
double fTolerance
 
bool fApplyFilter
 
int fMultMin [6]
 
int fMultMax [6]
 
TString fInvMassFilter
 
double fInvMassMin
 
double fInvMassMax
 
int fInvMassMult
 
int fCombIndex [5]
 
int fCombMult
 
bool fChargeConj
 
int fNAccept
 
bool fPersist
 
PndFsmDetFactoryfDetFac
 
std::string fAddedDets
 
FsmAbsDetList fDetList
 
TDatabasePDG * fdbPdg
 

Static Private Attributes

static TMatrixD fRho
 
static TMatrixD fEta
 

Detailed Description

Definition at line 32 of file PndFastSim.h.

Constructor & Destructor Documentation

PndFastSim::PndFastSim ( bool  persist = true)

Default constructor

Definition at line 62 of file PndFastSim.cxx.

References AddDetector(), fAddedDets, fApplyFilter, fChargeConj, fCombIndex, fCombMult, fdbPdg, fDetFac, fElectronBrems, fGenSplitOffs, fInvMassFilter, fInvMassMax, fInvMassMin, fInvMassMult, fMergeNeutralClusters, fMergeProbPar, fMultMax, fMultMin, fNAccept, fPersist, fPropagate, fRand, fToStartVtx, fUseCovMatrix, fUseFlatCovMatrix, fVb, i, and PndFsmRandom::Instance().

62  :
63 FairTask("Panda Fast Simulation") {
64  // fCandidates = new TClonesArray("TParticle");
67  fAddedDets=" ";
68  fVb=0;
69  fGenSplitOffs=false;
71  fMergeProbPar = 0.389;
72  fElectronBrems = false;
73  fUseFlatCovMatrix=true;
74  fPropagate=false;
75  fToStartVtx=false;
76  fUseCovMatrix=false;
77  fdbPdg = TDatabasePDG::Instance();
78  AddDetector("IdealPid");
79 
80  fInvMassFilter="";
81  fInvMassMin=0.;
82  fInvMassMax=1000.;
83  fInvMassMult=0;
84  fChargeConj=false;
85  fCombMult=0;
86  for (int i=0;i<5;++i) fCombIndex[i]=0;
87  for (int i=0;i<6;++i) {fMultMin[i]=0; fMultMax[i]=1000;}
88 
89 
90  fApplyFilter=false;
91  fNAccept = 0;
92 
93  fPersist = persist;
94 }
int fCombIndex[5]
Definition: PndFastSim.h:131
Int_t i
Definition: run_full.C:25
double fInvMassMin
Definition: PndFastSim.h:128
bool AddDetector(std::string name, std::string params="")
Definition: PndFastSim.cxx:313
bool fMergeNeutralClusters
Definition: PndFastSim.h:114
std::string fAddedDets
Definition: PndFastSim.h:139
bool fElectronBrems
Definition: PndFastSim.h:115
bool fPropagate
Definition: PndFastSim.h:117
int fMultMin[6]
Definition: PndFastSim.h:125
TString fInvMassFilter
Definition: PndFastSim.h:127
bool fApplyFilter
Definition: PndFastSim.h:124
bool fPersist
Definition: PndFastSim.h:136
TRandom3 * fRand
Definition: PndFastSim.h:108
int fNAccept
Definition: PndFastSim.h:134
double fMergeProbPar
Definition: PndFastSim.h:116
bool fGenSplitOffs
Definition: PndFastSim.h:113
static TRandom3 * Instance()
Definition: PndFsmRandom.cxx:4
bool fUseCovMatrix
Definition: PndFastSim.h:119
PndFsmDetFactory * fDetFac
Definition: PndFastSim.h:138
int fCombMult
Definition: PndFastSim.h:132
bool fUseFlatCovMatrix
Definition: PndFastSim.h:120
TDatabasePDG * fdbPdg
Definition: PndFastSim.h:142
double fInvMassMax
Definition: PndFastSim.h:129
int fInvMassMult
Definition: PndFastSim.h:130
bool fChargeConj
Definition: PndFastSim.h:133
int fMultMax[6]
Definition: PndFastSim.h:126
bool fToStartVtx
Definition: PndFastSim.h:118
PndFastSim::~PndFastSim ( )

Destructor

Definition at line 98 of file PndFastSim.cxx.

References fDetFac, fEventInfo, fMcCandidates, fPidChargedCand, fPidChargedProb, fPidNeutralCand, fPidNeutralProb, and fRand.

98  {
99  FairRootManager *fManager =FairRootManager::Instance();
100  fManager->Write();
101 
102  //if (fChargedCandidates) {fChargedCandidates->Delete(); delete fChargedCandidates;}
103  //if (fNeutralCandidates) {fNeutralCandidates->Delete(); delete fNeutralCandidates;}
104  if (fMcCandidates) {fMcCandidates->Delete(); delete fMcCandidates;}
105  //if (fMicroCandidates) {fMicroCandidates->Delete(); delete fMicroCandidates;}
106  if (fPidChargedCand) {fPidChargedCand->Delete(); delete fPidChargedCand;}
107  if (fPidNeutralCand) {fPidNeutralCand->Delete(); delete fPidNeutralCand;}
108 
109  if (fPidChargedProb) {fPidChargedProb->Delete(); delete fPidChargedProb;}
110  if (fPidNeutralProb) {fPidNeutralProb->Delete(); delete fPidNeutralProb;}
111 
112  if (fEventInfo) {fEventInfo->Delete(); delete fEventInfo;}
113 
114  if (fDetFac) delete fDetFac;
115  if (fRand) delete fRand;
116 }
TClonesArray * fMcCandidates
Definition: PndFastSim.h:97
TRandom3 * fRand
Definition: PndFastSim.h:108
TClonesArray * fPidNeutralProb
PndPidProbability TCA for charged particles.
Definition: PndFastSim.h:102
TClonesArray * fPidNeutralCand
Definition: PndFastSim.h:99
TClonesArray * fEventInfo
PndPidProbability TCA&#39;s for individual detectors.
Definition: PndFastSim.h:106
PndFsmDetFactory * fDetFac
Definition: PndFastSim.h:138
TClonesArray * fPidChargedCand
Definition: PndFastSim.h:98
TClonesArray * fPidChargedProb
Definition: PndFastSim.h:101

Member Function Documentation

int PndFastSim::acceptFilters ( RhoCandList l)
private

Definition at line 429 of file PndFastSim.cxx.

References RhoCandList::Add(), chCon(), RhoCandList::Combine(), RhoCandList::CombineAndAppend(), copyAndSetMass(), d, fabs(), fChargeConj, fCombIndex, fCombMult, fInvMassMax, fInvMassMin, fInvMassMult, fMultMax, fMultMin, fVb, RhoCandList::GetLength(), i, idx, and RhoCandList::Select().

Referenced by Exec().

430 {
431  RhoCandList plus, minus,lsp[14],comb, combsel;
432 
433  int i,j, nplus, nminus, ngam;
434 
435  for (i=0;i<l.GetLength();++i)
436  {
437  if (l[i]->Charge()>0) plus.Add(l[i]);
438  else if (l[i]->Charge()<0) minus.Add(l[i]);
439  else if (fabs(l[i]->Charge())<1e-6) lsp[10].Add(l[i]);
440  }
441 
442  // apply multiplicity filters
443 
444  // positiv particles
445  nplus = plus.GetLength();
446  if (nplus<fMultMin[0] || nplus>fMultMax[0]) return 1;
447 
448  // negative particles
449  nminus = minus.GetLength();
450  if (nminus<fMultMin[1] || nminus>fMultMax[1]) return 2;
451 
452  // neutrals
453  ngam = lsp[10].GetLength();
454  if (ngam<fMultMin[2] || ngam>fMultMax[2]) return 3;
455 
456  int pdgcodes[10] = {-11, 11, -13, 13, 211, -211, 321, -321, 2212, -2212};
457 
459  RhoMassParticleSelector pi0sel("pi0sel",0.135,0.06);
460  RhoMassParticleSelector kssel("kssel",0.4976,0.08);
461  RhoMassParticleSelector etasel("msel",0.547,0.08);
462 
463  // pi0s
464  if (fMultMin[3]>0 || fMultMax[3]<1000)
465  {
466  lsp[11].Combine(lsp[10],lsp[10]);
467  lsp[11].Select(&pi0sel);
468 
469  int npi0=lsp[11].GetLength();
470  if (npi0<fMultMin[3] || npi0>fMultMax[3]) return 4;
471  }
472 
473  // ks
474  if (fMultMin[4]>0 || fMultMax[4]<1000)
475  {
476 
477  copyAndSetMass(plus, lsp[4], pdgcodes[4]);
478  copyAndSetMass(minus, lsp[5], pdgcodes[5]);
479 
480  lsp[12].Combine(lsp[4],lsp[5]);
481  lsp[12].Select(&kssel);
482 
483  int nks=lsp[12].GetLength();
484  if (nks<fMultMin[4] || nks>fMultMax[4]) return 5;
485  }
486 
487  // etas
488  if (fMultMin[5]>0 || fMultMax[5]<1000)
489  {
490  lsp[13].Combine(lsp[10],lsp[10]);
491  lsp[13].Select(&etasel);
492 
493  int neta=lsp[13].GetLength();
494  if (neta<fMultMin[5] || neta>fMultMax[5]) return 6;
495  }
496 
497  // apply inv mass filter
498  // list indices: 0:e+ 1:e-, 2:mu+, ... 9:p-, 10:gamma, 11:pi0(gg), 12:KS, 13:eta(gg)
499  if (fInvMassMult>0)
500  {
501  if (fVb>0) cout <<"Comb mult = "<<fCombMult<<endl;
502  for (i=0;i<fCombMult;++i)
503  {
504  // prepare all needed lists
505  int idx = fCombIndex[i];
506  int ccidx = chCon(idx);
507 
508  if (fVb) cout <<"IDX:"<<idx<<" CCIDX:"<<ccidx<<endl;
509 
510  // charged particle species
511  if (idx<10 && lsp[idx].GetLength()==0)
512  {
513  if (idx%2) copyAndSetMass(minus, lsp[idx], pdgcodes[idx]);
514  else copyAndSetMass(plus, lsp[idx], pdgcodes[idx]);
515  }
516 
517  // if charged conjugation needed, prepare those lists
518  if (fChargeConj && ccidx<10 && lsp[ccidx].GetLength()==0)
519  {
520  if (ccidx%2) copyAndSetMass(minus, lsp[ccidx], pdgcodes[ccidx]);
521  else copyAndSetMass(plus, lsp[ccidx], pdgcodes[ccidx]);
522  }
523 
524  // pi0s
525  if (idx==11 && lsp[idx].GetLength()==0)
526  {
527  lsp[11].Combine(lsp[10],lsp[10]);
528  lsp[11].Select(&pi0sel);
529  }
530 
531  // ks
532  if (idx==12 && lsp[idx].GetLength()==0)
533  {
534  if (lsp[4].GetLength()==0) copyAndSetMass(plus, lsp[4], pdgcodes[4]);
535  if (lsp[5].GetLength()==0) copyAndSetMass(minus, lsp[5], pdgcodes[5]);
536 
537  lsp[12].Combine(lsp[4],lsp[5]);
538  lsp[12].Select(&kssel);
539  }
540 
541  // eta
542  if (idx==13 && lsp[idx].GetLength()==0)
543  {
544  lsp[13].Combine(lsp[10],lsp[10]);
545  lsp[13].Select(&etasel);
546  }
547  }
548 
549  int *iar=fCombIndex;
550 
551  switch (fCombMult)
552  {
553  case 2: comb.Combine(lsp[iar[0]], lsp[iar[1]]);
554  if (fChargeConj)
555  comb.CombineAndAppend(lsp[chCon(iar[0])], lsp[chCon(iar[1])]);
556  break;
557  case 3: comb.Combine(lsp[iar[0]], lsp[iar[1]], lsp[iar[2]]);
558  if (fChargeConj)
559  comb.CombineAndAppend(lsp[chCon(iar[0])], lsp[chCon(iar[1])], lsp[chCon(iar[2])]);
560  break;
561  case 4: comb.Combine(lsp[iar[0]], lsp[iar[1]], lsp[iar[2]], lsp[iar[3]]);
562  if (fChargeConj)
563  comb.CombineAndAppend(lsp[chCon(iar[0])], lsp[chCon(iar[1])], lsp[chCon(iar[2])], lsp[chCon(iar[3])]);
564  break;
565  case 5: comb.Combine(lsp[iar[0]], lsp[iar[1]], lsp[iar[2]], lsp[iar[3]], lsp[iar[4]]);
566  if (fChargeConj)
567  comb.CombineAndAppend(lsp[chCon(iar[0])], lsp[chCon(iar[1])], lsp[chCon(iar[2])],
568  lsp[chCon(iar[3])], lsp[chCon(iar[4])]);
569  break;
570  }
571 
572  combsel.Select(comb,&msel);
573  if (combsel.GetLength()<fInvMassMult)
574  {
575  if (fVb>0)
576  {
577  cout<<"filter list masses: ";
578  for (i=0;i<comb.GetLength();++i)
579  {
580  std::cout <<comb[i]->M()<<" ("<<comb[i]->GetMarker(0)<<") "<<endl;
581  cout <<" ("<<comb[i]->P4().X()<<","<<comb[i]->P4().Y()<<","<<comb[i]->P4().Z()<<","<<comb[i]->P4().E()<<")"<<endl;
582  for (j=0;j<comb[i]->NDaughters();++j)
583  {
584  RhoCandidate *d=comb[i]->Daughter(j);
585  cout<<*d<<endl;
586  }
587  }
588  cout<<endl;
589  }
590  return 7;
591  }
592 
593  }
594 
595  return 0;
596 }
int fCombIndex[5]
Definition: PndFastSim.h:131
int chCon(int i)
Definition: PndFastSim.cxx:422
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
TObjArray * d
Int_t i
Definition: run_full.C:25
double fInvMassMin
Definition: PndFastSim.h:128
Int_t GetLength() const
Definition: RhoCandList.h:46
int fMultMin[6]
Definition: PndFastSim.h:125
int idx[MAX]
Definition: autocutx.C:38
void CombineAndAppend(RhoCandList &l1, RhoCandList &l2)
void Combine(RhoCandList &l1, RhoCandList &l2)
void Select(RhoParticleSelectorBase *pidmgr)
void copyAndSetMass(RhoCandList &l, RhoCandList &nl, int pdg)
Definition: PndFastSim.cxx:362
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
int fCombMult
Definition: PndFastSim.h:132
double fInvMassMax
Definition: PndFastSim.h:129
int fInvMassMult
Definition: PndFastSim.h:130
bool fChargeConj
Definition: PndFastSim.h:133
int fMultMax[6]
Definition: PndFastSim.h:126
bool PndFastSim::AddDetector ( std::string  name,
std::string  params = "" 
)

Definition at line 313 of file PndFastSim.cxx.

References PndFsmDetFactory::create(), fAddedDets, fDetFac, fDetList, and fVb.

Referenced by PndFastSim(), prod_fsim(), QAmacro_fastsim_1(), quickfsimana(), run_fast(), simfast(), simfast_cmp(), simfast_dpm(), simfast_dpm_cmp(), simfast_opt(), simfast_single_allpid(), simfast_singletracks(), and tut_fastsim().

314 {
315  if (name =="")
316  {
317  cout <<" -W- (PndFastSim::AddDetector) - No name given, no detector added!"<<endl;
318 
319  }
320 
321  if (fAddedDets.find(" "+name+" ")!=std::string::npos)
322  {
323  std::cout <<" -W- (PndFastSim::AddDetector) - Detector <"<<name<<"> already appended - skipping!"<<std::endl;
324  return false;
325  }
326  fAddedDets.append(std::string(name+" "));
327 
328  PndFsmAbsDet *det=0;
329  det=fDetFac->create(name,params);
330  if (!det) return false;
331 
332  fDetList.push_back(det);
333  if (fVb>0) std::cout<<" -I- (PndFastSim::AddDetector) - Added detector "<<name<<" with params <"<<params<<">"<<std::endl;
334  return true;
335 }
std::string fAddedDets
Definition: PndFastSim.h:139
PndFsmAbsDet * create(std::string &name, ArgList &par)
TString name
FsmAbsDetList fDetList
Definition: PndFastSim.h:140
PndFsmDetFactory * fDetFac
Definition: PndFastSim.h:138
bool PndFastSim::AddDetector ( PndFsmAbsDet det)

Definition at line 337 of file PndFastSim.cxx.

References PndFsmAbsDet::detName(), fAddedDets, fDetList, and fVb.

338 {
339  if (det) {
340  fDetList.push_back(det);
341  fAddedDets.append(det->detName()+" ");
342  if (fVb>0) std::cout<<" -I- (PndFastSim::AddDetector) - Added detector "<<det->detName()<<std::endl;
343  }
344  return det;
345 }
std::string fAddedDets
Definition: PndFastSim.h:139
const std::string & detName()
Definition: PndFsmAbsDet.h:74
FsmAbsDetList fDetList
Definition: PndFastSim.h:140
int PndFastSim::chCon ( int  i)
private

Definition at line 422 of file PndFastSim.cxx.

References i.

Referenced by acceptFilters().

423 {
424  int cc=i;
425  if (i<10) i%2 ? cc=i-1 : cc=i+1;
426  return cc;
427 }
Int_t i
Definition: run_full.C:25
PndFastSim::ClassDef ( PndFastSim  ,
 
)
private
void PndFastSim::copyAndSetMass ( RhoCandList l,
RhoCandList nl,
int  pdg 
)
private

Definition at line 362 of file PndFastSim.cxx.

References RhoCandList::Add(), c, RhoCandList::Cleanup(), RhoCandList::GetLength(), i, and RhoCandidate::SetType().

Referenced by acceptFilters().

363 {
364  nl.Cleanup();
365  for (int i=0;i<l.GetLength();++i)
366  {
367  RhoCandidate c(*(l[i]));
368  c.SetType(pdg);
369  nl.Add(&c);
370  }
371 }
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
void Cleanup()
Definition: RhoCandList.cxx:62
Int_t i
Definition: run_full.C:25
Int_t GetLength() const
Definition: RhoCandList.h:46
bool PndFastSim::cutAndSmear ( PndFsmTrack t,
PndFsmResponse r 
)
private

Definition at line 1039 of file PndFastSim.cxx.

References c, PndFsmTrack::charge(), dE, PndFsmResponse::dE(), PndFsmResponse::detected(), PndFsmResponse::dm(), PndFsmResponse::dp(), PndFsmResponse::dphi(), dtheta, PndFsmResponse::dtheta(), PndFsmResponse::dV(), fabs(), fEta, fPropagate, fRand, fRho, fTolerance, fToStartVtx, fUseCovMatrix, fUseFlatCovMatrix, PndFsmTrack::GetHelixCov(), PndFsmTrack::GetHelixOmega(), PndFsmTrack::GetHelixParams(), PndFsmTrack::GetHelixTanDip(), PndFsmTrack::HelixRep(), ir, PndFsmResponse::m2(), m2(), PndFsmResponse::MvddEdx(), p, p2, PndFsmTrack::p4(), PndFsmTrack::Propagate(), SetFlatCovMatrix(), sin(), smearEnergy(), smearM(), smearM2(), smearMomentum(), smearMvddEdx(), smearPhi(), smearSttdEdx(), smearTheta(), smearVertex(), sqrt(), PndFsmTrack::startVtx(), PndFsmResponse::SttdEdx(), theta, and UpdateGammaHit().

Referenced by cutAndSmear(), and smearTrack().

1040 {
1041  if( !(r->detected()) ) return false;
1042 
1043  double charge = t->charge();
1044  double dE = r->dE();
1045  double dp = r->dp();
1046  double dtheta = r->dtheta();
1047  double dphi = r->dphi();
1048  double dm = r->dm();
1049  double m2 = r->m2();
1050  double MvddEdx =r->MvddEdx();
1051  double SttdEdx =r->SttdEdx();
1052  TVector3 dV = r->dV();
1053  //r->print(cout);
1054  //this removes candidates, which only have hit a PID-only device (like Cherenkov, TOF ...)
1055  if (fabs(charge)>1e-6 && fabs(dp)<1e-8) return false;
1056  if (fabs(charge)<1e-6 && fabs(dE)<1e-8) return false;
1057 
1058  // now produce some correlated errors
1059  if (fPropagate && fabs(charge)>1e-6) {
1060  t->HelixRep(t->startVtx() );
1061  double p2 = t->p4().Vect().Mag2();
1062  double omega = t->GetHelixOmega();
1063  double tandip = t->GetHelixTanDip();
1064  double theta = t->p4().Theta();
1065 
1066  static TVectorD gaus(5);
1067  for (char p=0;p<5;p++) gaus[p]=fRand->Gaus();
1068  if (fUseCovMatrix) gaus*=fEta;
1069  // calculate track par errors
1070  double err[5];
1071  // d0 (guessed), phi0, z0
1072  err[0] = 0.5*(dV.X()+dV.Y());
1073  err[1] = dphi;
1074  err[3] = dV.Z();
1075  // omega needs some error propagation (momentum and theta are uncorrelated)
1076  err[2] = omega * sqrt( dp*dp/p2 + pow(dtheta*tandip,2) );
1077  // as well as tandip
1078  err[4] = dtheta/pow(sin(theta),2);
1079  // smear track pars
1080  for (int p=0;p<5;p++) t->GetHelixParams()[p] += err[p] * gaus[p];
1081 
1082  // write scaled cov matrix
1083  for (int ir=0;ir<5;ir++) {
1084  for (int c=0;c<5;c++) {
1085  t->GetHelixCov()(ir,c) = fRho(ir,c)*fabs(err[ir]*err[c])*(fUseCovMatrix||ir==c);
1086  }
1087  }
1088  t->Propagate( fToStartVtx ? t->startVtx() : TVector3(0,0,0), fTolerance);
1089  // uncharged particles remain uncorrelated
1090  } else {
1091  //first set cov, using mctruth for correct calculation of jacobian
1092  if (fUseFlatCovMatrix) SetFlatCovMatrix(t,dp,dtheta,dphi,dE,dV.X(),dV.Y(),dV.Z());
1093  if (dE != 0.0) smearEnergy(t,dE);
1094  if (dp != 0.0) smearMomentum(t,dp);
1095  if (dtheta != 0.0) smearTheta(t,dtheta);
1096  if (dphi != 0.0) smearPhi(t,dphi);
1097  if ( fabs(charge)>1e-6 && (dV.X() != 0.0 || dV.Y() != 0.0 || dV.Z() != 0.0)) smearVertex(t,dV);
1098  if ( fabs(charge)<1e-6 ) UpdateGammaHit(t);
1099  }
1100  if (dm != 0.0) smearM(t,dm);
1101  if( m2!=0.0) smearM2(t,m2); // mass^2 of track after tof
1102  if(MvddEdx!=0.0) smearMvddEdx(t,MvddEdx); // dEdx of track after Mvd
1103  if(SttdEdx!=0.0) smearSttdEdx(t,SttdEdx); // dEdx of track after Stt
1104  return true;
1105 }
TVector3 dV()
void smearMomentum(PndFsmTrack *t, double dp)
void SetFlatCovMatrix(PndFsmTrack *t, double dp=0., double dtheta=0., double dphi=0., double dE=0., double dx=0., double dy=0., double dz=0.)
Double_t p
Definition: anasim.C:58
void smearPhi(PndFsmTrack *t, double dphi)
double fTolerance
Definition: PndFastSim.h:121
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t GetHelixOmega() const
Definition: PndFsmTrack.h:139
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
void smearSttdEdx(PndFsmTrack *t, double dedx)
double MvddEdx()
double charge()
Definition: PndFsmTrack.h:75
static TMatrixD fEta
Definition: PndFastSim.h:152
Double_t dE
Definition: anasim.C:58
TMatrixD & GetHelixCov()
Definition: PndFsmTrack.h:131
bool fPropagate
Definition: PndFastSim.h:117
double * GetHelixParams()
Definition: PndFsmTrack.h:130
TString m2(TString pts, TString exts="e px py pz")
Definition: invexp.C:117
void smearMvddEdx(PndFsmTrack *t, double dedx)
TLorentzVector p4()
Definition: PndFsmTrack.h:72
double SttdEdx()
void smearM2(PndFsmTrack *t, double m2)
double dtheta
Definition: anaLmdCluster.C:54
void smearTheta(PndFsmTrack *t, double dtheta)
void smearEnergy(PndFsmTrack *t, double dE)
TRandom3 * fRand
Definition: PndFastSim.h:108
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
void UpdateGammaHit(PndFsmTrack *t)
void smearM(PndFsmTrack *t, double dm)
TPad * p2
Definition: hist-t7.C:117
TVector3 startVtx()
Definition: PndFsmTrack.h:73
void HelixRep(TVector3 reference)
Definition: PndFsmTrack.cxx:98
bool fUseCovMatrix
Definition: PndFastSim.h:119
bool fUseFlatCovMatrix
Definition: PndFastSim.h:120
Double_t GetHelixTanDip() const
Definition: PndFsmTrack.h:141
void smearVertex(PndFsmTrack *t, TVector3 dV)
void Propagate(TVector3 origin, double deltaError=2.5)
static TMatrixD fRho
Definition: PndFastSim.h:151
bool fToStartVtx
Definition: PndFastSim.h:118
static int ir
Definition: ranlxd.cxx:374
bool PndFastSim::cutAndSmear ( PndFsmTrack t)
private

Definition at line 1031 of file PndFastSim.cxx.

References cutAndSmear(), and PndFsmTrack::detResponse().

1032 {
1033  return cutAndSmear(t, t->detResponse());
1034 }
bool cutAndSmear(PndFsmTrack *t, PndFsmResponse *r)
PndFsmResponse * detResponse()
Definition: PndFsmTrack.h:78
void PndFastSim::EnableElectronBremsstrahlung ( bool  brems = true)
inline

Definition at line 62 of file PndFastSim.h.

References fElectronBrems.

Referenced by prod_fsim(), quickfsimana(), simfast(), simfast_opt(), and tut_fastsim().

62 {fElectronBrems=brems;}
bool fElectronBrems
Definition: PndFastSim.h:115
void PndFastSim::EnablePropagation ( bool  propagate = true,
bool  tostartvtx = true,
bool  usecovmatrix = true,
double  tolerance = 0.0 
)

Definition at line 303 of file PndFastSim.cxx.

References fPropagate, fTolerance, fToStartVtx, fUseCovMatrix, and propagate().

Referenced by QAmacro_fastsim_1(), run_fast(), simfast_cmp(), and simfast_dpm_cmp().

303  {
305  fToStartVtx=tostartvtx;
306  fUseCovMatrix=usecovmatrix;
307  fTolerance=tolerance;
308 }
void propagate(TLorentzVector &l, TVector3 &p, float charge, TH2F *hpro=0)
double fTolerance
Definition: PndFastSim.h:121
bool fPropagate
Definition: PndFastSim.h:117
bool fUseCovMatrix
Definition: PndFastSim.h:119
bool fToStartVtx
Definition: PndFastSim.h:118
bool PndFastSim::EnableSplitoffs ( std::string  fname = "splitpars.dat")

Definition at line 224 of file PndFastSim.cxx.

References fGenSplitOffs, fspo, fVb, and i.

Referenced by prod_fsim(), quickfsimana(), simfast(), simfast_dpm(), simfast_opt(), simfast_single_allpid(), simfast_singletracks(), and tut_fastsim().

225 {
226  if (fname=="")
227  {
228  cout <<" -W- (PndFastSim::EnableSplitoffs) - no filename given; no split offs will be produced."<<endl;
229  return false;
230  }
231 
232  /*
233  prepare the split off parametrization
234  expected in the parameter file:
235  - four blocks (mom, phi, tht, multiplicity) with the parametrization structured like
236  * <n_i> <lowlimit> <uplimit> <model-string(in ROOT TF1 format)>
237  * 5 lines of parameter sets with info for pid=11,13,211,321,2212
238  <pid> <p1> <p2> ... <p_n_i>
239 
240  for the moment momentum independent
241  models:
242  momentum : expo(0)+expo(2) = 4 params
243  tht+phi : gaus(0)+gaus(3) = 6 params
244  multiplicity : expo(0)+gaus(2) = 5 params
245  */
246 
247  ifstream pars(fname.data());
248 
249  if ( (pars.rdstate() & ifstream::failbit ) != 0 )
250  {
251  cout << " -W- (PndFastSim::EnableSplitoffs) - Error opening '"<<fname<<"'; no split offs will be produced."<<endl;
252  return false;
253  }
254 
255 
256  int i,j,k;
257  char tmp[50];
258 
259  string model[4];
260  int numpars[4];
261  double rangeup[4],rangelow[4];
262  double curpar;
263 
264 
265  //loop over blocks
266  for (i=0;i<4;i++)
267  {
268  pars>>model[i]>>numpars[i]>>rangelow[i]>>rangeup[i];
269  if (fVb) cout <<" -I- (PndFastSim::EnableSplitoffs) split off model #"<<i<<": '"<<model[i]<<"' with "<<numpars[i]<<" params."<<endl;
270 
271  //loop over pid
272  for (j=0;j<5;j++)
273  {
274  sprintf(tmp,"f%d%d",j,i);
275  fspo[j][i]=new TF1(tmp,model[i].data());
276 
277  //loop over params
278  pars>>curpar; // read the dummy PID value
279  for (k=0;k<numpars[i];k++)
280  {
281  pars>>curpar;
282  fspo[j][i]->SetParameter(k,curpar);
283  }
284 
285  fspo[j][i]->SetRange(rangelow[i],rangeup[i]);
286  if (fVb) fspo[j][i]->Print();
287  }
288  }
289 
290  if (fVb) cout <<" -I- (PndFastSim::EnableSplitoffs) - Successfully read "<<fname<<endl;
291 
292  pars.close();
293 
294  fGenSplitOffs=true;
295 
296  return true;
297 }
Int_t i
Definition: run_full.C:25
TF1 * fspo[5][4]
Definition: PndFastSim.h:111
bool fGenSplitOffs
Definition: PndFastSim.h:113
void PndFastSim::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 616 of file PndFastSim.cxx.

References acceptFilters(), RhoCandList::Add(), c, PndFsmTrack::charge(), PndFsmTrack::Cov7(), PndFsmTrack::detResponse(), PndFsmResponse::DrcBarrelThtc(), PndFsmResponse::DrcBarrelThtcErr(), PndFsmResponse::DrcDiscThtc(), PndFsmResponse::DrcDiscThtcErr(), PndFsmResponse::EmcEcal(), evtcnt, exp(), fabs(), fApplyFilter, fBremsEnergy, fdbPdg, fElectronBrems, fEventInfo, fGenSplitOffs, fMcCandidates, fMergeNeutralClusters, fMergeProbPar, fNAccept, fPidArrayList, fPidChargedCand, fPidChargedProb, fPidNeutralCand, fPidNeutralProb, fPropagate, fRand, fspo, fUseFlatCovMatrix, fVb, RhoCandList::GetLength(), PndStack::GetNtrack(), PndStack::GetParticle(), i, RhoFactory::Instance(), PndFsmResponse::LHElectron(), PndFsmResponse::LHKaon(), PndFsmResponse::LHMuon(), PndFsmResponse::LHPion(), PndFsmResponse::LHProton(), PndFsmResponse::m2(), mom, PndFsmResponse::MuoIron(), PndFsmResponse::MvddEdx(), PndFsmTrack::p4(), phot, pos, PndStack::PushTrack(), RhoFactory::Reset(), PndFsmResponse::RichThtc(), PndFsmResponse::RichThtcErr(), PndEventInfo::SetCharged(), PndEventInfo::SetCmFrame(), PndPidCandidate::SetCov7(), PndPidCandidate::SetDiscNumberOfPhotons(), PndPidCandidate::SetDiscThetaC(), PndPidCandidate::SetDiscThetaCErr(), PndPidCandidate::SetDrcNumberOfPhotons(), PndPidCandidate::SetDrcThetaC(), PndPidCandidate::SetDrcThetaCErr(), PndPidProbability::SetElectronPdf(), PndPidCandidate::SetEmcCalEnergy(), PndPidCandidate::SetEnergy(), PndPidProbability::SetIndex(), PndEventInfo::SetIPTruth(), PndPidProbability::SetKaonPdf(), PndPidCandidate::SetLastHit(), RhoCandidate::SetMarker(), PndPidCandidate::SetMcIndex(), RhoCandidate::SetMcTruth(), PndPidCandidate::SetMomentum(), PndPidCandidate::SetMuoIron(), PndPidProbability::SetMuonPdf(), PndPidCandidate::SetMvdDEDX(), PndEventInfo::SetNeutrals(), PndPidProbability::SetPionPdf(), PndPidProbability::SetProtonPdf(), PndPidCandidate::SetRichNumberOfPhotons(), PndPidCandidate::SetRichThetaC(), PndPidCandidate::SetRichThetaCErr(), PndPidCandidate::SetSttMeanDEDX(), PndPidCandidate::SetTofM2(), RhoCandidate::SetTrackNumber(), RhoCandidate::SetType(), smearTrack(), PndFsmTrack::startVtx(), PndFsmTrack::stopVtx(), PndFsmResponse::SttdEdx(), and t.

617 {
618  int nCharged = 0;
619  int nNeutral = 0;
620 
621  if ((++evtcnt)%100==0)
622  cout <<"[PndFastSim] evt "<<evtcnt<<endl;
623 
624  PndStack *fStack=(PndStack*)gMC->GetStack();
625  Int_t nTracks=fStack->GetNtrack();
626 
628  // Reset output array
629  if (fMcCandidates->GetEntriesFast() != 0) fMcCandidates->Clear("C");
630  if (fPidChargedCand->GetEntriesFast() != 0) fPidChargedCand->Clear("C");
631  if (fPidNeutralCand->GetEntriesFast() != 0) fPidNeutralCand->Clear("C");
632  if (fPidChargedProb->GetEntriesFast() != 0) fPidChargedProb->Clear("C");
633  if (fPidNeutralProb->GetEntriesFast() != 0) fPidNeutralProb->Clear("C");
634  //if (fMicroCandidates->GetEntriesFast() != 0) fMicroCandidates->Clear("C");
635  if (fEventInfo->GetEntriesFast() != 0) fEventInfo->Clear("C");
636  for(std::map<TString,TClonesArray*>::iterator iter=fPidArrayList.begin();iter!=fPidArrayList.end();iter++)
637  {
638  if(iter->second->GetEntriesFast() != 0) iter->second->Clear("C");
639  }
640 
641  TClonesArray &mctracks = *fMcCandidates;
642 
643  TClonesArray &chrgCandidates = *fPidChargedCand; // Charged Candidates
644  TClonesArray &neutCandidates = *fPidNeutralCand; // Neutral Candidates
645 
646  TClonesArray &chrgProbs = *fPidChargedProb; // PID for charged Cands
647  TClonesArray &neutProbs = *fPidNeutralProb; // PID for neutral Cands
648 
649  //TClonesArray &microCandidates = *fMicroCandidates;
650  TClonesArray &evtInfo = *fEventInfo;
651 
652  if (fVb)cout <<"number of tracks **** "<< nTracks <<endl;
653 
654 
655  RhoCandList l, lfilt;
656  TLorentzVector McSumP4(0,0,0,0);
657  TVector3 McAvgVtx(0,0,0);
658 
659  if (fApplyFilter)
660  {
661  // extra loop for the filters
662  for (Int_t iTrack=0; iTrack<nTracks; iTrack++)
663  {
664  TParticle *t = fStack->GetParticle(iTrack);
665  int pdg = abs(t->GetPdgCode());
666  if (!(pdg==11 || pdg==13 || pdg==211 || pdg==321 || pdg==2212 || pdg==22)) continue;
667 
668  TLorentzVector p4(t->Px(),t->Py(),t->Pz(),t->Energy());
669  TParticlePDG* part = fdbPdg->GetParticle(t->GetPdgCode());
670  double charge=0;//safety against unknown pdgcodes, might scrw up the charge
671  if (part) charge=part->Charge();
672  if (fabs(charge)>2) charge/=3.;
673 
674  RhoCandidate c(p4, charge);
675  c.SetType(t->GetPdgCode());
676  c.SetMarker(lfilt.GetLength());
677  lfilt.Add(&c);
678  }
679 
680  int errac = acceptFilters(lfilt);
681 
682  if (errac>0)
683  {
684  if (fVb>0)
685  {
686  cout <<"PndFastSim Filter reject ev="<< evtcnt<<" with code="<<errac<<endl;
687  cout <<"Filter list"<<endl;
688  for (int i=0;i<lfilt.GetLength();++i)
689  {
690  cout <<i<<" : "<<lfilt[i]->PdgCode();
691  cout <<" ("<<lfilt[i]->P4().X()<<","<<lfilt[i]->P4().Y()<<","<<lfilt[i]->P4().Z()<<","<<lfilt[i]->P4().E()<<")"<<endl;
692  }
693  }
694 
696  return;
697  }
698  fNAccept++;
699  } // filter done
700 
701  for (Int_t iTrack=0; iTrack<nTracks; iTrack++)
702  {
703 
704  // Get the MCTrack information
705  Int_t mcsize = mctracks.GetEntriesFast();
706  Int_t chcandsize = chrgCandidates.GetEntriesFast();
707  Int_t neucandsize = neutCandidates.GetEntriesFast();
708  //Int_t miccandsize = microCandidates.GetEntriesFast();
709 
710  TParticle *t = fStack->GetParticle(iTrack);
711  if (fVb>1) t->Print();
712 
713  TLorentzVector p4(t->Px(),t->Py(),t->Pz(),t->Energy());
714  TVector3 stvtx(t->Vx(),t->Vy(),t->Vz());
715 
716  int pdg = t->GetPdgCode();
717 
718  // simulate bremsstrahlung for electrons and add photon on stack
719  if (abs(pdg)==11 && fElectronBrems)
720  {
721  // probability for bremsstrahlung was estimated from Full Sim to about 32%
722  if (fRand->Rndm()<0.32)
723  {
724  // get random energy loss and compute residual momentum (as equivalent to kinetic energy)
725  double loss = fBremsEnergy->GetRandom(0.03,0.9);
726 
727  // modify electron momentum mag (not the direction at the moment)
728  p4.SetVectM(p4.Vect()*(1.0-loss),5.11e-4);
729  TLorentzVector phot;
730  phot.SetVectM(p4.Vect()*loss, 0.0);
731 
732  // add an additional photon to the stack with the energy
733  fStack->PushTrack(0, iTrack, 22, // Int_t toBeDone, Int_t parentID, Int_t pdgCode
734  phot.X(), phot.Y(), phot.Z(), // Double_t px, Double_t py, Double_t pz,
735  phot.E(), 0., 0., // Double_t e, Double_t vx, Double_t vy,
736  0. , 0., 0., // Double_t vz, Double_t time, Double_t polx,
737  0., 0., kPPrimary, // Double_t poly, Double_t polz, TMCProcess proc,
738  nTracks, 0., 0.); // Int_t& ntr, Double_t weight, Int_t is
739 
740  nTracks=fStack->GetNtrack();
741  }
742  }
743 
744  // shall we add some parametrized split offs?
745  if (fGenSplitOffs)
746  {
747  int type=-1;
748 
749  if (abs(pdg) == 11) type=0;
750  else if (abs(pdg) == 211) type=2;
751  else if (abs(pdg) == 321) type=3;
752  else if (abs(pdg) == 2212) type=4;
753 
754  if (type>0)
755  {
756 
757  //number of split offs?
758  int numSP=(int)fspo[type][3]->GetRandom();
759 
760  if (fVb) cout <<" -I- (PndFastSim::Exec) - creating "<<numSP
761  <<" split offs for particle with type "<<type<<endl;
762 
763  for (int ii=0;ii<numSP;ii++)
764  {
765  TLorentzVector lv(p4);
766 
767  double mom = fspo[type][0]->GetRandom();
768  double dphi = fspo[type][1]->GetRandom();
769  double dtht = fspo[type][2]->GetRandom();
770 
771  lv.SetPhi(lv.Phi()+dphi);
772  lv.SetTheta(lv.Theta()+dtht);
773  lv.SetRho(mom);
774  lv.SetE(lv.P());
775 
776  // add an additional photon to the stack with the energy
777  fStack->PushTrack(0, iTrack, 22, // Int_t toBeDone, Int_t parentID, Int_t pdgCode
778  lv.X(), lv.Y(), lv.Z(), // Double_t px, Double_t py, Double_t pz,
779  lv.E(), 0., 0., // Double_t e, Double_t vx, Double_t vy,
780  0. , 0., 0., // Double_t vz, Double_t time, Double_t polx,
781  0., 0., kPPrimary, // Double_t poly, Double_t polz, TMCProcess proc,
782  nTracks, 0., 0.); // Int_t& ntr, Double_t weight, Int_t is
783 
784  nTracks=fStack->GetNtrack();
785  } // split off loop
786  }
787  }// generate split offs
788 
789  //TLorentzVector vtx(stvtx,t->T());
790  TParticlePDG* part = fdbPdg->GetParticle(t->GetPdgCode());
791  double charge(0);//safety against unknown pdgcodes, might scrw up the charge
792  if(part) charge=part->Charge();
793  if (fabs(charge)>2) charge/=3.;
794 
795  PndFsmTrack *ft=new PndFsmTrack(p4,stvtx,stvtx,charge,t->GetPdgCode(),iTrack+1);
796  //PndFsmTrack ft(p4,stvtx,stvtx,charge,t->GetPdgCode(),iPoint+1);
797 
798  // store a plain copy of the mc track to the file
799  RhoCandidate *pmc=new (mctracks[mcsize]) RhoCandidate(ft->p4(),ft->charge());
800  pmc->SetMcTruth(pmc);
801  pmc->SetPos(ft->startVtx());
802  pmc->SetType(t->GetPdgCode());
803  // write some ideal pid lhs
804  bool pdgcheck = false;
805  switch(abs(t->GetPdgCode())) {
806  case 22: pdgcheck = true; break;
807  case 11: pmc->SetPidInfo(0, 1); pdgcheck = true; break;
808  case 13: pmc->SetPidInfo(1, 1); pdgcheck = true; break;
809  case 211: pmc->SetPidInfo(2, 1); pdgcheck = true; break;
810  case 321: pmc->SetPidInfo(3, 1); pdgcheck = true; break;
811  case 2212: pmc->SetPidInfo(4, 1); pdgcheck = true; break;
812  }
813 
814  McSumP4+=ft->p4();
815  McAvgVtx+=ft->startVtx();
816 
817  if (pdgcheck) // only consider final state particles
818  {
819  // smear and cut the track according to the detector setup
820  bool smeared = smearTrack(ft, chcandsize);
821  if (smeared)
822  {
823  RhoVector3Err *svtx=new RhoVector3Err(ft->startVtx());
824  TLorentzVector miclv=ft->p4();
825  TVector3 pos=ft->startVtx();
826  TVector3 lastpos=ft->stopVtx();
827 
828  // assign pion mass to all charged and 0 to all neutral cands
829  if (fabs(ft->charge())>0.001)
830  {
831  miclv.SetVectM(miclv.Vect(),fdbPdg->GetParticle(211)->Mass());
832  nCharged++;
833  }
834  else
835  {
836  miclv.SetVectM(miclv.Vect(),0.0);
837  nNeutral++;
838  }
839 
840  PndPidCandidate *pidCand=0;
841  PndPidProbability *pidProb=0;
842 
843  if (fabs(ft->charge())>1e-6)
844  {
845  pidCand = new (chrgCandidates[chcandsize]) PndPidCandidate((Int_t)ft->charge(),pos,miclv);
846  pidProb = new (chrgProbs[chcandsize]) PndPidProbability();
847  }
848  else
849  {
850  // flag for merge neutral clusters
851  bool merged = false;
852 
853  // if two neutrals have small angle difference alpha, merge their clusters with a certain probability 1-exp(-a*alpha)
855  {
856  // loop through old gammas
857  for (int i=0; i<neucandsize;++i)
858  {
859  TLorentzVector nlv = ((PndPidCandidate*) neutCandidates[i])->GetLorentzVector();
860  double openang = nlv.Angle(miclv.Vect())*57.2958;
861  if (fRand->Rndm()>(1-exp(-fMergeProbPar*openang)))
862  {
863  double mergeE = nlv.E()+miclv.E();
864  TVector3 mergeV = nlv.Vect()+miclv.Vect();
865  mergeV *= mergeE/mergeV.Mag();
866 
867  PndPidCandidate* mergedCand = (PndPidCandidate*) neutCandidates[i];
868  mergedCand->SetMomentum(mergeV);
869  mergedCand->SetEnergy(mergeE);
870  mergedCand->SetMcIndex(-1); // remove MC truth match
871 
872  merged = true;
873 
874  // exit loop, if mergeing happend
875  i=neucandsize;
876  }
877  }
878  }
879 
880  // if the new gamma wasn't merged to another, just add it as usual
881  if (!merged)
882  {
883  pidCand = new (neutCandidates[neucandsize]) PndPidCandidate((Int_t)ft->charge(),pos,miclv);
884  pidCand->SetLastHit(lastpos);
885  pidProb = new (neutProbs[neucandsize]) PndPidProbability();
886  }
887  }
888 
889  if (pidProb)
890  {
891  pidProb->SetElectronPdf(ft->detResponse()->LHElectron());
892  pidProb->SetMuonPdf(ft->detResponse()->LHMuon());
893  pidProb->SetPionPdf(ft->detResponse()->LHPion());
894  pidProb->SetKaonPdf(ft->detResponse()->LHKaon());
895  pidProb->SetProtonPdf(ft->detResponse()->LHProton());
896  pidProb->SetIndex(chcandsize);
897  }
898  if (pidCand)
899  {
900  pidCand->SetMcIndex(iTrack);
901  pidCand->SetMvdDEDX( ft->detResponse()->MvddEdx() );
902 
903  //pidCand->SetMvdDEdxErr( ft->detResponse()->MvddEdxErr() );
904  pidCand->SetSttMeanDEDX( ft->detResponse()->SttdEdx() );
905  //pidCand->SetSttDEdxErr( ft->detResponse()->SttdEdxErr() );
906  pidCand->SetTofM2( ft->detResponse()->m2() );
907  //pidCand->SetTofM2Err( ft->detResponse()->m2Err() );
908  pidCand->SetDrcThetaC( ft->detResponse()->DrcBarrelThtc() );
909  pidCand->SetDrcThetaCErr( ft->detResponse()->DrcBarrelThtcErr() );
910  pidCand->SetDrcNumberOfPhotons(0);
911  pidCand->SetDiscThetaC( ft->detResponse()->DrcDiscThtc() );
912  pidCand->SetDiscThetaCErr( ft->detResponse()->DrcDiscThtcErr() );
913  pidCand->SetDiscNumberOfPhotons(0);
914  pidCand->SetRichThetaC( ft->detResponse()->RichThtc() );
915  pidCand->SetRichThetaCErr( ft->detResponse()->RichThtcErr() );
916  pidCand->SetRichNumberOfPhotons(0);
917  pidCand->SetEmcCalEnergy(ft->detResponse()->EmcEcal() );
918  pidCand->SetMuoIron(ft->detResponse()->MuoIron() );
919  RhoCandidate tcand(ft->p4(),ft->charge(),svtx);
920  tcand.SetTrackNumber(chcandsize);
921 
922  if( (fPropagate||fUseFlatCovMatrix))// && fabs(ft->charge())>1e-6)
923  {
924  tcand.SetCov7(ft->Cov7());
925  pidCand->SetCov7( ft->Cov7() );
926  }
927  l.Add(&tcand);
928  }
929  delete svtx;
930 
931  }// smeartrack
932  }
933  delete ft;
934  }//trackloop
935 
936  McAvgVtx*=1./(double)nTracks;
937 
938  PndEventInfo *eventInfo=new (evtInfo[evtInfo.GetEntriesFast()]) PndEventInfo();
939 
940  eventInfo->SetIPTruth(McAvgVtx);
941  eventInfo->SetCmFrame(McSumP4);
942  eventInfo->SetCharged(nCharged);
943  eventInfo->SetNeutrals(nNeutral);
944 
945  // TEventShape shape(l);
946  // eventInfo->SetEventShape(shape);
947 
948 }
TVector3 pos
void SetMomentum(TVector3 &mom)
double MuoIron()
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
void SetNeutrals(int n)
Definition: PndEventInfo.h:66
std::map< TString, TClonesArray * > fPidArrayList
PndPidProbability TCA for neutral particles.
Definition: PndFastSim.h:104
TLorentzVector phot
TVector3 stopVtx()
Definition: PndFsmTrack.h:74
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
TClonesArray * fMcCandidates
Definition: PndFastSim.h:97
Int_t i
Definition: run_full.C:25
void SetPionPdf(Double_t val)
int acceptFilters(RhoCandList &l)
Definition: PndFastSim.cxx:429
virtual Int_t GetNtrack() const
Definition: PndStack.h:116
Int_t GetLength() const
Definition: RhoCandList.h:46
static void Reset()
Definition: RhoFactory.cxx:28
void SetEmcCalEnergy(Double_t val)
void SetCov7(const TMatrixD &cov7)
bool fMergeNeutralClusters
Definition: PndFastSim.h:114
TF1 * fspo[5][4]
Definition: PndFastSim.h:111
void SetMcIndex(int idx)
TParticle * GetParticle(Int_t trackId) const
Definition: PndStack.cxx:442
void SetKaonPdf(Double_t val)
double MvddEdx()
double charge()
Definition: PndFsmTrack.h:75
double RichThtcErr()
void SetTofM2(Double_t val)
Double_t mom
Definition: plot_dirc.C:14
void SetDiscThetaC(Double_t val)
void SetCmFrame(TLorentzVector &cmf)
double RichThtc()
double DrcDiscThtcErr()
void SetElectronPdf(Double_t val)
void SetCharged(int n)
Definition: PndEventInfo.h:65
void SetMuonPdf(Double_t val)
bool fElectronBrems
Definition: PndFastSim.h:115
bool fPropagate
Definition: PndFastSim.h:117
void SetEnergy(Double_t en)
void SetDrcThetaCErr(Double_t val)
TMatrixD & Cov7()
Definition: PndFsmTrack.h:132
TLorentzVector p4()
Definition: PndFsmTrack.h:72
double SttdEdx()
bool fApplyFilter
Definition: PndFastSim.h:124
void SetMuoIron(Double_t val)
void SetDrcNumberOfPhotons(Int_t val)
void SetRichNumberOfPhotons(Int_t val)
double LHElectron()
TRandom3 * fRand
Definition: PndFastSim.h:108
int fNAccept
Definition: PndFastSim.h:134
double fMergeProbPar
Definition: PndFastSim.h:116
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
static RhoFactory * Instance()
Definition: RhoFactory.cxx:34
bool fGenSplitOffs
Definition: PndFastSim.h:113
TClonesArray * fPidNeutralProb
PndPidProbability TCA for charged particles.
Definition: PndFastSim.h:102
void SetSttMeanDEDX(Double_t val)
void SetProtonPdf(Double_t val)
bool smearTrack(PndFsmTrack *t, int idx=-1)
Definition: PndFastSim.cxx:951
void SetIPTruth(TVector3 &inVtx)
Definition: PndEventInfo.h:58
void SetMcTruth(RhoCandidate *mct)
Definition: RhoCandidate.h:436
double DrcDiscThtc()
TClonesArray * fPidNeutralCand
Definition: PndFastSim.h:99
TVector3 startVtx()
Definition: PndFsmTrack.h:73
void SetTrackNumber(Int_t trnum=-1)
Definition: RhoCandidate.h:418
void SetRichThetaC(Double_t val)
void SetDrcThetaC(Double_t val)
void SetDiscNumberOfPhotons(Int_t val)
TClonesArray * fEventInfo
PndPidProbability TCA&#39;s for individual detectors.
Definition: PndFastSim.h:106
TF1 * fBremsEnergy
Definition: PndFastSim.h:112
void SetLastHit(TVector3 &pos)
TTree * t
Definition: bump_analys.C:13
bool fUseFlatCovMatrix
Definition: PndFastSim.h:120
virtual void PushTrack(Int_t toBeDone, Int_t parentID, Int_t pdgCode, Double_t px, Double_t py, Double_t pz, Double_t e, Double_t vx, Double_t vy, Double_t vz, Double_t time, Double_t polx, Double_t poly, Double_t polz, TMCProcess proc, Int_t &ntr, Double_t weight, Int_t is)
Definition: PndStack.cxx:57
TDatabasePDG * fdbPdg
Definition: PndFastSim.h:142
TClonesArray * fPidChargedCand
Definition: PndFastSim.h:98
PndFsmResponse * detResponse()
Definition: PndFsmTrack.h:78
TClonesArray * fPidChargedProb
Definition: PndFastSim.h:101
double DrcBarrelThtcErr()
void SetDiscThetaCErr(Double_t val)
double EmcEcal()
double DrcBarrelThtc()
void SetMvdDEDX(Double_t val)
void SetIndex(Int_t idx)
void SetRichThetaCErr(Double_t val)
void PndFastSim::Finish ( )
virtual

Definition at line 599 of file PndFastSim.cxx.

References evtcnt, fApplyFilter, fInvMassFilter, fInvMassMax, fInvMassMin, fInvMassMult, fMultMax, fMultMin, fNAccept, i, and TString.

600 {
601  TString multnames[6] = {"ch+","ch-", "gam", "pi0", "ks", "eta"};
602  if (fApplyFilter)
603  {
604  for (int i=0;i<6;++i)
605  {
606  if (fMultMin[i]>0 || fMultMax[i]<1000)
607  std::cout <<"*** Filtering Info: "<<fMultMin[i]<<" <= "<<multnames[i].Data() <<" <= "<<fMultMax[i]<<endl;
608  }
609  std::cout <<"*** Filtering Info: "<<fInvMassMin<<" < m("<<fInvMassFilter.Data()<<") < "<<fInvMassMax<<" with N="<<fInvMassMult<<endl;
610  std::cout <<"*** Filtering Info: "<< fNAccept <<"/"<<evtcnt<<" accepted by filters."<<endl;
611  }
612 }
Int_t i
Definition: run_full.C:25
double fInvMassMin
Definition: PndFastSim.h:128
int fMultMin[6]
Definition: PndFastSim.h:125
TString fInvMassFilter
Definition: PndFastSim.h:127
bool fApplyFilter
Definition: PndFastSim.h:124
int fNAccept
Definition: PndFastSim.h:134
double fInvMassMax
Definition: PndFastSim.h:129
int fInvMassMult
Definition: PndFastSim.h:130
int fMultMax[6]
Definition: PndFastSim.h:126
InitStatus PndFastSim::Init ( )
virtual

Virtual method Init

Definition at line 158 of file PndFastSim.cxx.

References evtcnt, fBremsEnergy, fDetList, fEta, fRho, fVb, PndFsmAbsDet::print(), and Register().

158  {
159 
160  if (fVb>3) cout << " Inside the Init function****" << endl;
161 
162  Register();
163  if (fVb)
164  {
165  cout <<"\n\n****** DETECTOR SETUP SUMMARY ***************\n\n";
166  for (FsmAbsDetList::iterator iter=fDetList.begin();iter!=fDetList.end(); iter++)
167  {
168  cout<<endl<<"----------------------------"<<endl;
169  PndFsmAbsDet *det=*iter;
170  det->print(std::cout);
171  }
172  cout <<"\n\n****** DETECTOR SETUP SUMMARY ***************\n\n";
173  }
174 
175  evtcnt=0;
176 
177  // for the time being, this is the matrix of correlation
178  // coefficients, it is later scaled by the track errors
179  //d0 phi0 omega z0 dandip
180  fRho(0,0)= 1; fRho(0,1)=-0.9243; fRho(0,2)= 0.3087; fRho(0,3)=-0.01743; fRho(0,4)= 0.009165; // d0
181  fRho(1,0)=-0.9243; fRho(1,1)= 1; fRho(1,2)=-0.4189; fRho(1,3)=-0.00552; fRho(1,4)=-0.00873; // phi0
182  fRho(2,0)= 0.3087; fRho(2,1)=-0.4189; fRho(2,2)= 1; fRho(2,3)= 0.02094; fRho(2,4)=-0.01742; // omega
183  fRho(3,0)=-0.01743; fRho(3,1)=-0.00552; fRho(3,2)= 0.02094; fRho(3,3)= 1; fRho(3,4)=-0.9327; // z0
184  fRho(4,0)= 0.009165; fRho(4,1)=-0.00873; fRho(4,2)=-0.01742; fRho(4,3)=-0.9327; fRho(4,4)= 1; // tandip
185  // and this is the square root of fRho
186  // needed to correlate the track params
187  //d0 phi0 omega z0 dandip
188  fEta(0,0)= 1; // d0
189  fEta(1,0)=-0.9243; fEta(1,1)= 0.3817; // phi0
190  fEta(2,0)= 0.3087; fEta(2,1)=-0.35; fEta(2,2)= 0.8844; // omega
191  fEta(3,0)=-0.01743; fEta(3,1)=-0.05667; fEta(3,2)= 0.007335; fEta(3,3)= 0.9982; // z0
192  fEta(4,0)= 0.009165; fEta(4,1)=-6.781e-4; fEta(4,2)=-0.02316; fEta(4,3)=-0.9341; fEta(4,4)=0.3562; // tandip
193 
194  // parametrization for bremsstrahlung (based on Crystal Ball function tail)
195  fBremsEnergy = new TF1("fBremsEnergy","[0]*pow([3]/[4],[3])*exp(-0.5*[4]^2)/pow((x-[1])/[2]+[3]/[4]-[4],[3])");
196  fBremsEnergy->SetParameters(1,0.005,0.0311,1.05,1.04);
197  fBremsEnergy->SetRange(0.03,0.9); // valid between 3% and 90% energy loss through bremsstrahlung
198 
199  // Create and register output array
200  cout << "-I- PndFastSim: Intialization successfull" << endl;
201 
202  return kSUCCESS;
203 }
static TMatrixD fEta
Definition: PndFastSim.h:152
virtual void print(std::ostream &o)
FsmAbsDetList fDetList
Definition: PndFastSim.h:140
void Register()
Definition: PndFastSim.cxx:120
TF1 * fBremsEnergy
Definition: PndFastSim.h:112
static TMatrixD fRho
Definition: PndFastSim.h:151
bool PndFastSim::MergeNeutralClusters ( bool  merge = true,
double  par = 0.389 
)
inline

Definition at line 61 of file PndFastSim.h.

References fMergeNeutralClusters, fMergeProbPar, merge(), and par.

Referenced by prod_fsim(), quickfsimana(), simfast(), simfast_opt(), and tut_fastsim().

Double_t par[3]
bool fMergeNeutralClusters
Definition: PndFastSim.h:114
int merge(TString ntp="", TString fout="", TString f1="", TString f2="", TString f3="", TString f4="", TString f5="")
double fMergeProbPar
Definition: PndFastSim.h:116
void PndFastSim::Register ( )

Definition at line 120 of file PndFastSim.cxx.

References detname, PndFsmAbsDet::detName(), PndFsmAbsDet::doesPid(), fDetList, fEventInfo, fMcCandidates, fPersist, fPidArrayList, fPidChargedCand, fPidChargedProb, fPidNeutralCand, fPidNeutralProb, and TString.

Referenced by Init().

120  {
121  //---
122 
123  fMcCandidates = new TClonesArray("RhoCandidate");
124  FairRootManager::Instance()->Register("PndMcTracks","FastSim", fMcCandidates, fPersist);
125 
126  fPidChargedCand = new TClonesArray("PndPidCandidate");
127  FairRootManager::Instance()->Register("PidChargedCand","FastSim", fPidChargedCand, fPersist);
128 
129  fPidNeutralCand = new TClonesArray("PndPidCandidate");
130  FairRootManager::Instance()->Register("PidNeutralCand","FastSim", fPidNeutralCand, fPersist);
131 
132  fPidChargedProb = new TClonesArray("PndPidProbability");
133  FairRootManager::Instance()->Register("PidChargedProbability","FastSim", fPidChargedProb, fPersist);
134 
135  fPidNeutralProb = new TClonesArray("PndPidProbability");
136  FairRootManager::Instance()->Register("PidNeutralProbability","FastSim", fPidNeutralProb, fPersist);
137 
138  fEventInfo = new TClonesArray("PndEventInfo");
139  FairRootManager::Instance()->Register("PndEventSummary","FastSim", fEventInfo, fPersist);
140 
141  for (FsmAbsDetList::iterator iter=fDetList.begin();iter!=fDetList.end(); iter++)
142  {
143  PndFsmAbsDet *det=*iter;
144  if(!det->doesPid()) continue;
145  TString detname = det->detName();
146  TString arrayname=detname;
147  arrayname+="Probability";
148  TClonesArray* tmparray = new TClonesArray("PndPidProbability");
149  FairRootManager::Instance()->Register(arrayname.Data(),"FastSim", tmparray, fPersist);
150  fPidArrayList[detname]=tmparray;
151  //std::cout<<"registered pid TCA with name "<<detname.Data()<<std::endl;
152  }
153 }
std::map< TString, TClonesArray * > fPidArrayList
PndPidProbability TCA for neutral particles.
Definition: PndFastSim.h:104
TClonesArray * fMcCandidates
Definition: PndFastSim.h:97
TString detname
Definition: anasim.C:61
const std::string & detName()
Definition: PndFsmAbsDet.h:74
bool fPersist
Definition: PndFastSim.h:136
FsmAbsDetList fDetList
Definition: PndFastSim.h:140
TClonesArray * fPidNeutralProb
PndPidProbability TCA for charged particles.
Definition: PndFastSim.h:102
TClonesArray * fPidNeutralCand
Definition: PndFastSim.h:99
Bool_t doesPid() const
Definition: PndFsmAbsDet.h:76
TClonesArray * fEventInfo
PndPidProbability TCA&#39;s for individual detectors.
Definition: PndFastSim.h:106
TClonesArray * fPidChargedCand
Definition: PndFastSim.h:98
TClonesArray * fPidChargedProb
Definition: PndFastSim.h:101
void PndFastSim::SetFlatCovMatrix ( PndFsmTrack t,
double  dp = 0.,
double  dtheta = 0.,
double  dphi = 0.,
double  dE = 0.,
double  dx = 0.,
double  dy = 0.,
double  dz = 0. 
)
private

Definition at line 1109 of file PndFastSim.cxx.

References PndFsmTrack::charge(), cos(), dE, dtheta, dx, dy, dz, fabs(), p, PndFsmTrack::p4(), PndFsmTrack::SetP4Cov(), PndFsmTrack::SetVCov(), and sin().

Referenced by cutAndSmear().

1110 {
1111  TLorentzVector lv=t->p4();
1112 
1113  double st=sin(lv.Theta());
1114  double ct=cos(lv.Theta());
1115  double sf=sin(lv.Phi());
1116  double cf=cos(lv.Phi());
1117  double p=lv.P();
1118  double e=lv.E();
1119 
1120  //printf("FastSim covariance test: dp=%f, dtheta=%f, dphi=%f, dE=%f, dx=%f, dy=%f, dz=%f\n", dp, dtheta, dphi, dE, dx, dy, dz);
1121 
1122  TMatrixD jacobian(4,3);
1123 
1124  if(fabs(t->charge())>1e-6)
1125  {
1126  jacobian(0,0) = st*cf;
1127  jacobian(1,0) = st*sf;
1128  jacobian(2,0) = ct;
1129  jacobian(3,0) = (e>0) ? p/e : 0.;
1130  } else { // no direct momentum measurement of neutrals
1131  jacobian(0,0) = st*cf*e/p;
1132  jacobian(1,0) = st*sf*e/p;
1133  jacobian(2,0) = ct*e/p;
1134  jacobian(3,0) = 1.;
1135  }
1136  jacobian(0,1) = p*ct*cf;
1137  jacobian(1,1) = p*ct*sf;
1138  jacobian(2,1) = -p*st;
1139  jacobian(3,1) = 0.;
1140  jacobian(0,2) = -p*st*sf;
1141  jacobian(1,2) = p*st*cf;
1142  jacobian(2,2) = 0.;
1143  jacobian(3,2) = 0.;
1144 
1145  TMatrixDSym covPol(3);
1146  if(fabs(t->charge())>1e-6) covPol(0,0)=dp*dp;
1147  else covPol(0,0)=dE*dE;
1148 
1149  covPol(1,1)=dtheta*dtheta;
1150  covPol(2,2)=dphi*dphi;
1151 
1152  TMatrixD jcov(jacobian,TMatrixD::kMult,covPol);
1153  TMatrixD covCar(jcov,TMatrixD::kMultTranspose,jacobian);
1154 
1155 // cout<<"Lorentzvector: ";
1156 // lv.Print();
1157 // cout<<"Polar covariance: ";
1158 // covPol.Print();
1159 // cout<<"Jacobian: ";
1160 // jacobian.Print();
1161 // //cout<<"Jacobian*covPol: ";
1162 // //jcov.Print();
1163 // cout<<"Kartesian covariance: ";
1164 // covCar.Print();
1165 
1166  t->SetP4Cov(covCar);
1167 
1168  TMatrixD covV(3,3);
1169  covV(0,0) = dx*dx;
1170  covV(0,1) = 0.;
1171  covV(0,2) = 0.;
1172  covV(1,0) = 0.;
1173  covV(1,1) = dy*dy;
1174  covV(1,2) = 0.;
1175  covV(2,0) = 0.;
1176  covV(2,1) = 0.;
1177  covV(2,2) = dz*dz;
1178 
1179  t->SetVCov(covV);
1180 
1181  return;
1182 }
Double_t p
Definition: anasim.C:58
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double dy
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
double charge()
Definition: PndFsmTrack.h:75
Double_t dE
Definition: anasim.C:58
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void SetVCov(TMatrixD &vcov)
double dtheta
Definition: anaLmdCluster.C:54
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double dx
void SetP4Cov(TMatrixD &p4cov)
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
void PndFastSim::SetInvMassFilter ( TString  filter,
double  min,
double  max,
int  mult = 1 
)

Definition at line 373 of file PndFastSim.cxx.

References fApplyFilter, fChargeConj, fCombIndex, fCombMult, fInvMassFilter, fInvMassMax, fInvMassMin, fInvMassMult, i, idx, max(), min(), mult, and TString.

Referenced by prod_fsim(), quickfsimana(), and simfast_opt().

374 {
375  fInvMassFilter=filter;
376  fInvMassMin=min;
377  fInvMassMax=max;
379  fApplyFilter=true;
380  fCombMult = 0;
381 
382  TString wk=filter;
383  if (wk.EndsWith("cc"))
384  {
385  fChargeConj = true;
386  wk=wk(0,wk.Length()-3);
387  }
388  wk+=" ";
389 
390  // chop the configuration string in pieces
391  // coding is: 0+=e+, 1+=mu+, ... , 4+=p+ (& cc), gm=gamma
392  // list indices: 0:e+ 1:e-, 2:mu+, ... 9:p-, 10:gamma, 11:pi0(gg), 12:KS, 13:eta(gg)
393  while (wk.Length()>0 && fCombMult<5)
394  {
395  int delim = wk.Index(" ");
396  int idx=0;
397  TString tok=wk(0,delim);
398  cout <<tok<<" ";
399  if (tok=="gam") idx=10;
400  else if (tok=="pi0") idx=11;
401  else if (tok=="ks") idx=12;
402  else if (tok=="eta") idx=13;
403  else if (tok.BeginsWith("e")) idx=0;
404  else if (tok.BeginsWith("mu")) idx=2;
405  else if (tok.BeginsWith("pi")) idx=4;
406  else if (tok.BeginsWith("k")) idx=6;
407  else if (tok.BeginsWith("p")) idx=8;
408 
409  if (tok.EndsWith("-")) idx++;
410 
412  fCombMult++;
413 
414  wk=wk(delim+1,1000);
415  }
416  cout <<"Inv Mass Filter combines:";
417  for (int i=0;i<fCombMult;++i) cout <<fCombIndex[i]<<" ";
418  if (fChargeConj) cout << " with c.c.";
419  cout <<endl;
420 }
int fCombIndex[5]
Definition: PndFastSim.h:131
Int_t i
Definition: run_full.C:25
double fInvMassMin
Definition: PndFastSim.h:128
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
int idx[MAX]
Definition: autocutx.C:38
TString fInvMassFilter
Definition: PndFastSim.h:127
bool fApplyFilter
Definition: PndFastSim.h:124
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
int fCombMult
Definition: PndFastSim.h:132
double fInvMassMax
Definition: PndFastSim.h:129
int fInvMassMult
Definition: PndFastSim.h:130
bool fChargeConj
Definition: PndFastSim.h:133
Double_t mult
void PndFastSim::SetMultFilter ( TString  type,
int  min,
int  max = 1000 
)

Definition at line 347 of file PndFastSim.cxx.

References fMultMax, fMultMin, i, max(), and min().

Referenced by simfast_opt().

348 {
349  int i=0;
350  if (type=="+") i=0;
351  else if (type=="-") i=1;
352  else if (type=="gam") i=2;
353  else if (type=="pi0") i=3;
354  else if (type=="ks") i=4;
355  else if (type=="eta") i=5;
356  else cout <<" WARNING : Unknown type '"<<type.Data()<< "' for muliplicity filter!"<<endl;
357 
358  fMultMin[i] = min;
359  fMultMax[i] = max;
360 }
Int_t i
Definition: run_full.C:25
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
int fMultMin[6]
Definition: PndFastSim.h:125
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
int fMultMax[6]
Definition: PndFastSim.h:126
void PndFastSim::SetParContainers ( )
privatevirtual

Get parameter containers

Definition at line 205 of file PndFastSim.cxx.

References run.

205  {
206 
207  // Get run and runtime database
208  FairRun* run = FairRun::Instance();
209  if ( ! run ) Fatal("SetParContainers", "No analysis run");
210 
211  FairRuntimeDb* db = run->GetRuntimeDb();
212  if ( ! db ) Fatal("SetParContainers", "No runtime database");
213 
214 
215 }
Int_t run
Definition: autocutx.C:47
void PndFastSim::SetSeed ( unsigned int  seed = 65539)

Definition at line 218 of file PndFastSim.cxx.

References fRand.

218  {
219  fRand->SetSeed(seed);
220 }
unsigned int seed
TRandom3 * fRand
Definition: PndFastSim.h:108
void PndFastSim::SetUseFlatCov ( bool  v = true)
inline

Definition at line 64 of file PndFastSim.h.

References fUseFlatCovMatrix, and v.

Referenced by prod_fsim(), quickfsimana(), simfast(), simfast_opt(), and tut_fastsim().

__m128 v
Definition: P4_F32vec4.h:4
bool fUseFlatCovMatrix
Definition: PndFastSim.h:120
void PndFastSim::SetVerbosity ( int  vb)
inline
void PndFastSim::smearEnergy ( PndFsmTrack t,
double  dE 
)
private

Definition at line 1185 of file PndFastSim.cxx.

References PndFsmTrack::charge(), fabs(), fRand, m, PndFsmTrack::p4(), PndFsmTrack::setP4(), and sqrt().

Referenced by cutAndSmear().

1186 {
1187  TLorentzVector p4=t->p4();
1188  double m=t->p4().M();
1189 
1190  // if (t->pdt()->lundId()==22) //gammas are always on mass shell
1191  if (fabs(t->charge()) < 0.1 ) //gammas are always on mass shell
1192  {
1193  // double rnd = pRand->Gaus(0.,dE);
1194  double rnd = fRand->Gaus(0. , dE);
1195 
1196  double newE = rnd + p4.E();
1197 
1198  p4.SetVectM(p4.Vect()*(newE/p4.E()) , 0.0);
1199  }
1200  else
1201  {
1202  double newE = fRand->Gaus(0.0,dE) + p4.E();
1203  double newP = sqrt(newE*newE - m*m);
1204 
1205  p4.SetE(newE);
1206  p4.SetVectMag(p4.Vect(),newP);
1207  }
1208 
1209  t->setP4(p4);
1210 }
__m128 m
Definition: P4_F32vec4.h:28
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
double charge()
Definition: PndFsmTrack.h:75
Double_t dE
Definition: anasim.C:58
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void setP4(TLorentzVector l)
TRandom3 * fRand
Definition: PndFastSim.h:108
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
void PndFastSim::smearM ( PndFsmTrack t,
double  dm 
)
private

Definition at line 1242 of file PndFastSim.cxx.

References PndFsmTrack::p4(), and PndFsmTrack::setP4().

Referenced by cutAndSmear().

1243 {
1244  TLorentzVector p4=t->p4();
1245  // double newM = p4.m() + RandGauss::shoot(0.,dm);
1246  double newM=p4.M() + dm;
1247 
1248  p4.SetVectM(p4.Vect(),newM);
1249 
1250  t->setP4(p4);
1251 }
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void setP4(TLorentzVector l)
void PndFastSim::smearM2 ( PndFsmTrack t,
double  m2 
)
private

Definition at line 1254 of file PndFastSim.cxx.

References PndFsmTrack::setMass2().

Referenced by cutAndSmear().

1255 {
1256 
1257  t->setMass2(m2);
1258 }
void setMass2(double c)
TString m2(TString pts, TString exts="e px py pz")
Definition: invexp.C:117
void PndFastSim::smearMomentum ( PndFsmTrack t,
double  dp 
)
private

Definition at line 1213 of file PndFastSim.cxx.

References fRand, PndFsmTrack::p4(), and PndFsmTrack::setP4().

Referenced by cutAndSmear().

1214 {
1215  TLorentzVector p4=t->p4();
1216  double newP = p4.Vect().Mag() + fRand->Gaus(0.0,dp);
1217  p4.SetVectM(p4.Vect().Unit()*newP,t->p4().M());
1218  t->setP4(p4);
1219 }
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void setP4(TLorentzVector l)
TRandom3 * fRand
Definition: PndFastSim.h:108
void PndFastSim::smearMvddEdx ( PndFsmTrack t,
double  dedx 
)
private

Definition at line 1261 of file PndFastSim.cxx.

References PndFsmTrack::setMvddEdX().

Referenced by cutAndSmear().

1262 {
1263 
1264  t->setMvddEdX(MvddEdx);
1265 }
void setMvddEdX(double c)
void PndFastSim::smearPhi ( PndFsmTrack t,
double  dphi 
)
private

Definition at line 1232 of file PndFastSim.cxx.

References fRand, PndFsmTrack::p4(), and PndFsmTrack::setP4().

Referenced by cutAndSmear().

1233 {
1234  TLorentzVector p4=t->p4();
1235  double newPhi = p4.Phi() + fRand->Gaus(0.,dphi);
1236  p4.SetPhi(newPhi);
1237 
1238  t->setP4(p4);
1239 }
TLorentzVector p4()
Definition: PndFsmTrack.h:72
void setP4(TLorentzVector l)
TRandom3 * fRand
Definition: PndFastSim.h:108
void PndFastSim::smearSttdEdx ( PndFsmTrack t,
double  dedx 
)
private

Definition at line 1268 of file PndFastSim.cxx.

References PndFsmTrack::setSttdEdX().

Referenced by cutAndSmear().

1269 {
1270 
1271  t->setSttdEdX(SttdEdx);
1272 }
void setSttdEdX(double c)
void PndFastSim::smearTheta ( PndFsmTrack t,
double  dtheta 
)
private

Definition at line 1222 of file PndFastSim.cxx.

References fRand, PndFsmTrack::p4(), and PndFsmTrack::setP4().

Referenced by cutAndSmear().

1223 {
1224  TLorentzVector p4=t->p4();
1225  double newTheta = p4.Theta() + fRand->Gaus(0.,dtheta);
1226  p4.SetTheta(newTheta);
1227 
1228  t->setP4(p4);
1229 }
TLorentzVector p4()
Definition: PndFsmTrack.h:72
double dtheta
Definition: anaLmdCluster.C:54
void setP4(TLorentzVector l)
TRandom3 * fRand
Definition: PndFastSim.h:108
void PndFastSim::smearTpcdEdx ( PndFsmTrack t,
double  dedx 
)
private
bool PndFastSim::smearTrack ( PndFsmTrack t,
int  idx = -1 
)
private

Definition at line 951 of file PndFastSim.cxx.

References PndFsmTrack::charge(), cutAndSmear(), PndFsmResponse::detector(), detname, PndFsmAbsDet::detName(), PndFsmAbsDet::doesPid(), fabs(), fDetList, fPidArrayList, PndFsmResponse::LHElectron(), PndFsmResponse::LHKaon(), PndFsmResponse::LHMuon(), PndFsmResponse::LHPion(), PndFsmResponse::LHProton(), PndFsmAbsDet::respond(), PndFsmTrack::setDetResponse(), PndPidProbability::SetIndex(), mrfdata_8b_error::success, sumResponse(), and TString.

Referenced by Exec().

952 {
953 
954 
955  FsmResponseList responseList;
956 
957  for (FsmAbsDetList::iterator iter=fDetList.begin();iter!=fDetList.end(); iter++)
958  {
959  PndFsmAbsDet *det=*iter;
960  //if (!det) {cout <<"--------------------------> outch"<<endl;continue;}
961 
962  PndFsmResponse *respo=det->respond(t);
963 
964  if (respo) responseList.push_back(respo);
965  }
966 
967  PndFsmResponse *allResponse=sumResponse(responseList);
968 
969  t->setDetResponse(allResponse);
970 
971  bool success = cutAndSmear(t);
972 
973  if(success && fabs(t->charge())>1e-6)
974  { // create dummy PID info for all detectors, then fill it with reason later
975  for(std::map<TString, TClonesArray*>::iterator pidit=fPidArrayList.begin() ; pidit != fPidArrayList.end() ; pidit++)
976  {
977  TClonesArray* myPidarray = pidit->second;
978  if (!myPidarray) Error("PndFastSim::smearTrack","Missing PidProb Array: \"%s\"",pidit->first.Data());
979  if (myPidarray->GetEntriesFast() != chcandsize ) Warning("PndFastSim::smearTrack","unequal array sizes: cand array:%i prob array \"%s\" :%i",chcandsize,pidit->first.Data(),myPidarray->GetEntriesFast());
980  PndPidProbability *apidProb=new((*myPidarray)[chcandsize]) PndPidProbability();
981  apidProb->SetIndex(chcandsize);
982  }
983  }
984 
985  for (FsmResponseList::iterator riter=responseList.begin(); riter!=responseList.end();riter++)
986  {
987  PndFsmResponse *resp=*riter;
988  if (!resp) continue;
989 
990  if(success && fabs(t->charge())>1e-6&&resp->detector()->doesPid())
991  { //save PID information only if particle is stored
992  TString detname = resp->detector()->detName();
993  TClonesArray* myPidarray = fPidArrayList[detname];
994  if (!myPidarray) Error("PndFastSim::smearTrack","Failed accessing PidProb Array: \"%s\"",detname.Data());
995  PndPidProbability *pidProb=(PndPidProbability*)myPidarray->At(chcandsize);
996  if (!pidProb) Error("PndFastSim::smearTrack","Failed accessing PidProb Object number %i from array \"%s\"",chcandsize,detname.Data());
997  {
998  double rawLHe = resp->LHElectron();
999  double rawLHmu = resp->LHMuon();
1000  double rawLHpi = resp->LHPion();
1001  double rawLHK = resp->LHKaon();
1002  double rawLHp = resp->LHProton();
1003 
1004  double sumRaw = rawLHe+rawLHmu+rawLHpi+rawLHK+rawLHp;
1005 
1006  if (sumRaw!=0.)
1007  {
1008  pidProb->SetElectronPdf(rawLHe/sumRaw);
1009  pidProb->SetMuonPdf(rawLHmu/sumRaw);
1010  pidProb->SetPionPdf(rawLHpi/sumRaw);
1011  pidProb->SetKaonPdf(rawLHK/sumRaw);
1012  pidProb->SetProtonPdf(rawLHp/sumRaw);
1013  } else {
1014  pidProb->SetElectronPdf(0.2);
1015  pidProb->SetMuonPdf(0.2);
1016  pidProb->SetPionPdf(0.2);
1017  pidProb->SetKaonPdf(0.2);
1018  pidProb->SetProtonPdf(0.2);
1019  }
1020  }
1021  }
1022  // and now clean up!
1023  delete resp;
1024  }
1025  responseList.clear();
1026 
1027  return success;
1028 }
virtual PndFsmResponse * respond(PndFsmTrack *t)=0
std::map< TString, TClonesArray * > fPidArrayList
PndPidProbability TCA for neutral particles.
Definition: PndFastSim.h:104
void setDetResponse(PndFsmResponse *resp)
TString detname
Definition: anasim.C:61
double charge()
Definition: PndFsmTrack.h:75
PndFsmAbsDet * detector()
const std::string & detName()
Definition: PndFsmAbsDet.h:74
PndFsmResponse * sumResponse(FsmResponseList respList)
double LHElectron()
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
FsmAbsDetList fDetList
Definition: PndFastSim.h:140
Bool_t doesPid() const
Definition: PndFsmAbsDet.h:76
bool cutAndSmear(PndFsmTrack *t, PndFsmResponse *r)
std::list< PndFsmResponse * > FsmResponseList
Definition: PndFastSim.h:30
void SetIndex(Int_t idx)
void PndFastSim::smearVertex ( PndFsmTrack t,
TVector3  dV 
)
private

Definition at line 1275 of file PndFastSim.cxx.

References fRand, PndFsmTrack::setStartVtx(), and PndFsmTrack::startVtx().

Referenced by cutAndSmear().

1276 {
1277  TVector3 vtx = t->startVtx();
1278 
1279  double newX = vtx.X() + fRand->Gaus( 0. , dV.X() );
1280  double newY = vtx.Y() + fRand->Gaus( 0. , dV.Y() );
1281  double newZ = vtx.Z() + fRand->Gaus( 0. , dV.Z() );
1282 
1283  vtx.SetX(newX);
1284  vtx.SetY(newY);
1285  vtx.SetZ(newZ);
1286 
1287  t->setStartVtx(vtx);
1288 }
TRandom3 * fRand
Definition: PndFastSim.h:108
void setStartVtx(TVector3 v)
TVector3 startVtx()
Definition: PndFsmTrack.h:73
PndFsmResponse * PndFastSim::sumResponse ( FsmResponseList  respList)
private

Definition at line 1294 of file PndFastSim.cxx.

References dE, PndFsmResponse::dE(), PndFsmResponse::detected(), PndFsmResponse::detector(), PndFsmAbsDet::detName(), PndFsmResponse::dm(), PndFsmResponse::dp(), PndFsmResponse::dphi(), PndFsmResponse::DrcBarrelThtc(), PndFsmResponse::DrcBarrelThtcErr(), PndFsmResponse::DrcDiscThtc(), PndFsmResponse::DrcDiscThtcErr(), PndFsmResponse::dt(), dtheta, PndFsmResponse::dtheta(), PndFsmResponse::dV(), PndFsmResponse::EmcEcal(), fabs(), fVb, PndFsmResponse::LHElectron(), PndFsmResponse::LHKaon(), PndFsmResponse::LHMuon(), PndFsmResponse::LHPion(), PndFsmResponse::LHProton(), PndFsmResponse::m2(), m2(), PndFsmResponse::m2Err(), PndFsmResponse::MuoIron(), PndFsmResponse::MvddEdx(), PndFsmResponse::MvddEdxErr(), PndFsmResponse::print(), PndFsmResponse::RichThtc(), PndFsmResponse::RichThtcErr(), PndFsmResponse::setdE(), PndFsmResponse::setDetected(), PndFsmResponse::setdm(), PndFsmResponse::setdp(), PndFsmResponse::setdphi(), PndFsmResponse::setDrcBarrelThtc(), PndFsmResponse::setDrcDiscThtc(), PndFsmResponse::setdt(), PndFsmResponse::setdtheta(), PndFsmResponse::setdV(), PndFsmResponse::setEmcEcal(), PndFsmResponse::setLHElectron(), PndFsmResponse::setLHKaon(), PndFsmResponse::setLHMuon(), PndFsmResponse::setLHPion(), PndFsmResponse::setLHProton(), PndFsmResponse::setm2(), PndFsmResponse::setMuoIron(), PndFsmResponse::setMvddEdx(), PndFsmResponse::setRichThtc(), PndFsmResponse::setSttdEdx(), sqrt(), PndFsmResponse::SttdEdx(), PndFsmResponse::SttdEdxErr(), and val.

Referenced by smearTrack().

1295 {
1296 
1297  bool detected=false;
1298 
1299  double dE=0.0;
1300  double dp=0.0;
1301  double dtheta=0.0;
1302  double dphi=0.0;
1303  double dt=0.0;
1304  double dm=0.0;
1305 
1306  double m2=0;
1307  double MvddEdx=0;
1308  double SttdEdx=0;
1309  double DrcDiscThtc=0;
1310  double DrcBarrelThtc=0;
1311  double RichThtc=0;
1312  double EmcEcal=0;
1313  double MuoIron=0;
1314 
1315  double m2Err=0;
1316  double MvddEdxErr=0;
1317  double SttdEdxErr=0;
1318  double DrcDiscThtcErr=0;
1319  double DrcBarrelThtcErr=0;
1320  double RichThtcErr=0;
1321 
1322  double dVx=0.0;
1323  double dVy=0.0;
1324  double dVz=0.0;
1325 
1326  double LH_e=1.0;
1327  double LH_mu=1.0;
1328  double LH_pi=1.0;
1329  double LH_K=1.0;
1330  double LH_p=1.0;
1331 
1332  double val=0.0;
1333 
1334  PndFsmResponse *allResponse=new PndFsmResponse();
1335 
1336  if (fVb>3) cout <<" *** SINGLE RESPONSES ***"<<endl<<"--------------------------------------"<<endl;
1337 
1338  for (FsmResponseList::iterator riter=respList.begin(); riter!=respList.end();riter++)
1339  {
1340  PndFsmResponse *resp=*riter;
1341 
1342  if ( fVb>3 ) resp->print(cout);
1343 
1344  if (resp->detector()->detName()=="IdealPid") continue; // don't use ideal PID
1345 
1346  detected = detected | resp->detected();
1347 
1348  if (resp->detected())
1349  {
1350  if (fabs(val = resp->dE()) > 1e-8) dE += 1/(val*val);
1351  if (fabs(val = resp->dp()) > 1e-8) dp += 1/(val*val);
1352  if (fabs(val = resp->dtheta())> 1e-8) dtheta += 1/(val*val);
1353  if (fabs(val = resp->dphi()) > 1e-8) dphi += 1/(val*val);
1354 
1355  if (fabs(val = resp->dt()) > 1e-8) dt += val*val;
1356  if (fabs(val = resp->dm()) > 1e-8) dm =val;
1357  if (fabs (val = resp->m2()) > 1e-11) m2=val;
1358  if (fabs (val = resp->MvddEdx()) > 1e-11) MvddEdx=val;
1359  if (fabs (val = resp->SttdEdx()) > 1e-11) SttdEdx=val;
1360  if (fabs (val = resp->DrcDiscThtc()) > 1e-11) DrcDiscThtc=val;
1361  if (fabs (val = resp->DrcBarrelThtc()) > 1e-11) DrcBarrelThtc=val;
1362  if (fabs (val = resp->RichThtc()) > 1e-11) RichThtc=val;
1363  if (fabs (val = resp->EmcEcal()) > 1e-11) EmcEcal=val;
1364  if (fabs (val = resp->MuoIron()) > 1e-11) MuoIron=val;
1365 
1366  if (fabs (val = resp->m2Err()) > 1e-11) m2Err=val;
1367  if (fabs (val = resp->MvddEdxErr()) > 1e-11) MvddEdxErr=val;
1368  if (fabs (val = resp->SttdEdxErr()) > 1e-11) SttdEdxErr=val;
1369  if (fabs (val = resp->DrcDiscThtcErr()) > 1e-11) DrcDiscThtcErr=val;
1370  if (fabs (val = resp->DrcBarrelThtcErr()) > 1e-11) DrcBarrelThtcErr=val;
1371  if (fabs (val = resp->RichThtcErr()) > 1e-11) RichThtcErr=val;
1372 
1373  if (fabs (val = resp->dV().X()) > 1e-11) dVx += 1/(val*val);
1374  if (fabs (val = resp->dV().Y()) > 1e-11) dVy += 1/(val*val);
1375  if (fabs (val = resp->dV().Z()) > 1e-11) dVz += 1/(val*val);
1376 
1377  double rawLHe = resp->LHElectron();
1378  double rawLHmu = resp->LHMuon();
1379  double rawLHpi = resp->LHPion();
1380  double rawLHK = resp->LHKaon();
1381  double rawLHp = resp->LHProton();
1382 
1383  double sumRaw = rawLHe+rawLHmu+rawLHpi+rawLHK+rawLHp;
1384 
1385  if (sumRaw>0) {
1386  rawLHe /= sumRaw;
1387  rawLHmu /= sumRaw;
1388  rawLHpi /= sumRaw;
1389  rawLHK /= sumRaw;
1390  rawLHp /= sumRaw;
1391  LH_e *= rawLHe;
1392  LH_mu *= rawLHmu;
1393  LH_pi *= rawLHpi;
1394  LH_K *= rawLHK;
1395  LH_p *= rawLHp;
1396  } else {
1397  LH_e *= 0.2;
1398  LH_mu *= 0.2;
1399  LH_pi *= 0.2;
1400  LH_K *= 0.2;
1401  LH_p *= 0.2;
1402  }
1403 
1404  //here a weighted Likelihood evaluation has to be done
1405 
1406  /*
1407  if (val = rawLHe) LH_e == 0.0 ? LH_e = val : LH_e *= val;
1408  if (val = rawLHmu) LH_mu == 0.0 ? LH_mu = val : LH_mu *= val;
1409  if (val = rawLHpi) LH_pi == 0.0 ? LH_pi = val : LH_pi *= val;
1410  if (val = rawLHK) LH_K == 0.0 ? LH_K = val : LH_K *= val;
1411  if (val = rawLHp) LH_p == 0.0 ? LH_p = val : LH_p *= val;
1412  */
1413  }
1414  }
1415 
1416  double sumLH = LH_e + LH_mu + LH_pi + LH_K + LH_p;
1417 
1418  if (sumLH>0) {
1419  LH_e /= sumLH;
1420  LH_mu /= sumLH;
1421  LH_pi /= sumLH;
1422  LH_K /= sumLH;
1423  LH_p /= sumLH;
1424  } else {
1425  LH_e = 0.2;
1426  LH_mu = 0.2;
1427  LH_pi = 0.2;
1428  LH_K = 0.2;
1429  LH_p = 0.2;
1430  }
1431 
1432  allResponse->setDetected(detected);
1433  allResponse->setdE( dE>0. ? 1/sqrt(dE) : 0.0 );
1434  allResponse->setdp( dp>0. ? 1/sqrt(dp) : 0.0 );
1435  allResponse->setdtheta( dtheta>0. ? 1/sqrt(dtheta) : 0.0 );
1436  allResponse->setdphi( dphi>0. ? 1/sqrt(dphi) : 0.0 );
1437  allResponse->setdt( sqrt(dt) );
1438  allResponse->setdm(dm);
1439 
1440  allResponse->setm2(m2, m2Err);
1441  allResponse->setMvddEdx(MvddEdx,MvddEdxErr);
1442  allResponse->setSttdEdx(SttdEdx,SttdEdxErr);
1443 
1444  allResponse->setDrcDiscThtc(DrcDiscThtc,DrcDiscThtcErr);
1445  allResponse->setDrcBarrelThtc(DrcBarrelThtc,DrcBarrelThtcErr);
1446  allResponse->setRichThtc(RichThtc,RichThtcErr);
1447  allResponse->setEmcEcal(EmcEcal);
1448  allResponse->setMuoIron(MuoIron);
1449 
1450  if (dVx > 0.) dVx=1./sqrt(dVx); else dVx = 0.0;
1451  if (dVy > 0.) dVy=1./sqrt(dVy); else dVy = 0.0;
1452  if (dVz > 0.) dVz=1./sqrt(dVz); else dVz = 0.0;
1453 
1454  allResponse->setdV( dVx , dVy , dVz );
1455 
1456  allResponse->setLHElectron(LH_e);
1457  allResponse->setLHMuon(LH_mu);
1458  allResponse->setLHPion(LH_pi);
1459  allResponse->setLHKaon(LH_K);
1460  allResponse->setLHProton(LH_p);
1461 
1462 
1463  if (fVb>2) cout <<"--------------------------------------"<<endl<<"*** OVERALL RESPONSE ***"<<endl<<"--------------------------------------" <<endl;
1464  if (fVb>2) allResponse->print(cout);
1465  if (fVb>2) cout <<"--------------------------------------"<<endl;
1466 
1467  return allResponse;
1468 }
TVector3 dV()
double MuoIron()
void setDrcBarrelThtc(double val, double err=0)
void setdV(TVector3 v)
void setdphi(double val)
void setLHElectron(double val)
void setLHProton(double val)
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
double MvddEdxErr()
void setLHMuon(double val)
void print(std::ostream &o)
double MvddEdx()
double RichThtcErr()
PndFsmAbsDet * detector()
void setm2(double val, double err=0)
double RichThtc()
double DrcDiscThtcErr()
Double_t dE
Definition: anasim.C:58
void setMuoIron(double val)
TString m2(TString pts, TString exts="e px py pz")
Definition: invexp.C:117
void setdt(double val)
double SttdEdxErr()
double SttdEdx()
void setDrcDiscThtc(double val, double err=0)
const std::string & detName()
Definition: PndFsmAbsDet.h:74
double dtheta
Definition: anaLmdCluster.C:54
void setdp(double val)
double LHElectron()
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
void setEmcEcal(double val)
void setdE(double val)
double DrcDiscThtc()
void setRichThtc(double val, double err=0)
void setDetected(bool isdet)
void setdtheta(double val)
void setSttdEdx(double val, double err=0)
void setdm(double val)
double DrcBarrelThtcErr()
double EmcEcal()
void setMvddEdx(double val, double err=0)
double DrcBarrelThtc()
void setLHPion(double val)
void PndFastSim::UpdateGammaHit ( PndFsmTrack t)
private

Definition at line 1470 of file PndFastSim.cxx.

References fabs(), PndFsmTrack::p4(), PndFsmTrack::setStopVtx(), PndFsmTrack::stopVtx(), x, y, and z.

Referenced by cutAndSmear().

1471 {
1472  TVector3 hitpos=t->stopVtx();
1473  double z=hitpos.z();
1474  if(z==0.)
1475  { //barrel - TODO convention to set z=0 initially.
1476  if(fabs(t->p4().Perp()) < 1e-9) return;
1477  double perpscale = hitpos.Perp()/t->p4().Perp();
1478  double x = t->p4().Px()*perpscale;
1479  double y = t->p4().Py()*perpscale;
1480  z = t->p4().Pz()*perpscale;
1481  hitpos.SetXYZ(x,y,z);
1482  } else
1483  { // z-fixed disks
1484  double zscale = z/t->p4().Pz();
1485  double x = t->p4().Px()*zscale;
1486  double y = t->p4().Py()*zscale;
1487  hitpos.SetXYZ(x,y,z);
1488  }
1489  t->setStopVtx(hitpos);
1490 
1491 }
TVector3 stopVtx()
Definition: PndFsmTrack.h:74
void setStopVtx(TVector3 v)
TLorentzVector p4()
Definition: PndFsmTrack.h:72
Double_t z
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Double_t x
Double_t y

Member Data Documentation

int PndFastSim::evtcnt
private

Definition at line 110 of file PndFastSim.h.

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

std::string PndFastSim::fAddedDets
private

Definition at line 139 of file PndFastSim.h.

Referenced by AddDetector(), and PndFastSim().

bool PndFastSim::fApplyFilter
private

Definition at line 124 of file PndFastSim.h.

Referenced by Exec(), Finish(), PndFastSim(), and SetInvMassFilter().

TF1* PndFastSim::fBremsEnergy
private

Definition at line 112 of file PndFastSim.h.

Referenced by Exec(), and Init().

bool PndFastSim::fChargeConj
private

Definition at line 133 of file PndFastSim.h.

Referenced by acceptFilters(), PndFastSim(), and SetInvMassFilter().

int PndFastSim::fCombIndex[5]
private

Definition at line 131 of file PndFastSim.h.

Referenced by acceptFilters(), PndFastSim(), and SetInvMassFilter().

int PndFastSim::fCombMult
private

Definition at line 132 of file PndFastSim.h.

Referenced by acceptFilters(), PndFastSim(), and SetInvMassFilter().

TDatabasePDG* PndFastSim::fdbPdg
private

Definition at line 142 of file PndFastSim.h.

Referenced by Exec(), and PndFastSim().

PndFsmDetFactory* PndFastSim::fDetFac
private

Definition at line 138 of file PndFastSim.h.

Referenced by AddDetector(), PndFastSim(), and ~PndFastSim().

FsmAbsDetList PndFastSim::fDetList
private

Definition at line 140 of file PndFastSim.h.

Referenced by AddDetector(), Init(), Register(), and smearTrack().

bool PndFastSim::fElectronBrems
private

Definition at line 115 of file PndFastSim.h.

Referenced by EnableElectronBremsstrahlung(), Exec(), and PndFastSim().

TMatrixD PndFastSim::fEta
staticprivate

Definition at line 152 of file PndFastSim.h.

Referenced by cutAndSmear(), and Init().

TClonesArray* PndFastSim::fEventInfo
private

PndPidProbability TCA's for individual detectors.

Definition at line 106 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

bool PndFastSim::fGenSplitOffs
private

Definition at line 113 of file PndFastSim.h.

Referenced by EnableSplitoffs(), Exec(), and PndFastSim().

TString PndFastSim::fInvMassFilter
private

Definition at line 127 of file PndFastSim.h.

Referenced by Finish(), PndFastSim(), and SetInvMassFilter().

double PndFastSim::fInvMassMax
private

Definition at line 129 of file PndFastSim.h.

Referenced by acceptFilters(), Finish(), PndFastSim(), and SetInvMassFilter().

double PndFastSim::fInvMassMin
private

Definition at line 128 of file PndFastSim.h.

Referenced by acceptFilters(), Finish(), PndFastSim(), and SetInvMassFilter().

int PndFastSim::fInvMassMult
private

Definition at line 130 of file PndFastSim.h.

Referenced by acceptFilters(), Finish(), PndFastSim(), and SetInvMassFilter().

TClonesArray* PndFastSim::fMcCandidates
private

Output array of Candidates

Definition at line 97 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

bool PndFastSim::fMergeNeutralClusters
private

Definition at line 114 of file PndFastSim.h.

Referenced by Exec(), MergeNeutralClusters(), and PndFastSim().

double PndFastSim::fMergeProbPar
private

Definition at line 116 of file PndFastSim.h.

Referenced by Exec(), MergeNeutralClusters(), and PndFastSim().

TClonesArray* PndFastSim::fMicroCandidates
private

Definition at line 100 of file PndFastSim.h.

int PndFastSim::fMultMax[6]
private

Definition at line 126 of file PndFastSim.h.

Referenced by acceptFilters(), Finish(), PndFastSim(), and SetMultFilter().

int PndFastSim::fMultMin[6]
private

Definition at line 125 of file PndFastSim.h.

Referenced by acceptFilters(), Finish(), PndFastSim(), and SetMultFilter().

int PndFastSim::fNAccept
private

Definition at line 134 of file PndFastSim.h.

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

bool PndFastSim::fPersist
private

Definition at line 136 of file PndFastSim.h.

Referenced by PndFastSim(), and Register().

std::map<TString,TClonesArray*> PndFastSim::fPidArrayList
private

PndPidProbability TCA for neutral particles.

Definition at line 104 of file PndFastSim.h.

Referenced by Exec(), Register(), and smearTrack().

TClonesArray* PndFastSim::fPidChargedCand
private

Definition at line 98 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

TClonesArray* PndFastSim::fPidChargedProb
private

Definition at line 101 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

TClonesArray* PndFastSim::fPidNeutralCand
private

Definition at line 99 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

TClonesArray* PndFastSim::fPidNeutralProb
private

PndPidProbability TCA for charged particles.

Definition at line 102 of file PndFastSim.h.

Referenced by Exec(), Register(), and ~PndFastSim().

bool PndFastSim::fPropagate
private

Definition at line 117 of file PndFastSim.h.

Referenced by cutAndSmear(), EnablePropagation(), Exec(), and PndFastSim().

TRandom3* PndFastSim::fRand
private
TMatrixD PndFastSim::fRho
staticprivate

Definition at line 151 of file PndFastSim.h.

Referenced by cutAndSmear(), and Init().

TF1* PndFastSim::fspo[5][4]
private

Definition at line 111 of file PndFastSim.h.

Referenced by EnableSplitoffs(), and Exec().

double PndFastSim::fTolerance
private

Definition at line 121 of file PndFastSim.h.

Referenced by cutAndSmear(), and EnablePropagation().

bool PndFastSim::fToStartVtx
private

Definition at line 118 of file PndFastSim.h.

Referenced by cutAndSmear(), EnablePropagation(), and PndFastSim().

bool PndFastSim::fUseCovMatrix
private

Definition at line 119 of file PndFastSim.h.

Referenced by cutAndSmear(), EnablePropagation(), and PndFastSim().

bool PndFastSim::fUseFlatCovMatrix
private

Definition at line 120 of file PndFastSim.h.

Referenced by cutAndSmear(), Exec(), PndFastSim(), and SetUseFlatCov().

int PndFastSim::fVb
private

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