FairRoot/PandaRoot
PndPidDrcAssociatorTask.cxx
Go to the documentation of this file.
2 #include "PndPidCandidate.h"
3 #include "PndPidProbability.h"
4 #include "FairRootManager.h"
5 #include "TMath.h"
6 #include "TF1.h"
7 #include "Riostream.h"
8 
9 
10 //___________________________________________________________
12  //
13  FairRootManager *fManager =FairRootManager::Instance();
14  fManager->Write();
15 }
16 
17 //___________________________________________________________
19  //---
20  fPidChargedProb = new TClonesArray("PndPidProbability");
22 }
23 
24 //___________________________________________________________
25 PndPidDrcAssociatorTask::PndPidDrcAssociatorTask(const char *name, const char *title):FairTask(name)
26 {
27  //---
28  fPidChargedProb = new TClonesArray("PndPidProbability");
30  SetTitle(title);
31  }
32 
33 //___________________________________________________________
35 
36  std::cout << "InitStatus PndPidDrcAssociatorTask::Init()" << std::endl;
37 
38  FairRootManager *fManager =FairRootManager::Instance();
39 
40  fPidChargedCand = (TClonesArray *)fManager->GetObject("PidChargedCand"+fTrackBranchNamePidHypo);
41  if ( ! fPidChargedCand) {
42  std::cout << "-I- PndPidDrcAssociatorTask::Init: No PndPidCandidate array PidChargedCand there!" << std::endl;
43  return kERROR;
44  }
45 
46  Register();
47 
48  std::cout << "-I- PndPidDrcAssociatorTask::Init: Success!" << std::endl;
49 
50  return kSUCCESS;
51 }
52 
53 //______________________________________________________
55  //--
56 }
57 //______________________________________________________
59  if (fPidChargedProb->GetEntriesFast() != 0) fPidChargedProb->Clear();
60  if(fVerbose>1) std::cout << "-I- Start PndPidDrcAssociatorTask. "<<std::endl;
61 
62  // Get the Candidates
63  for(Int_t i=0; i<fPidChargedCand->GetEntriesFast(); i++)
64  {
66  TClonesArray& pidRef = *fPidChargedProb;
67  PndPidProbability* prob = new(pidRef[i]) PndPidProbability();// initializes with zeros
68  prob->SetIndex(i);
69  if (pidcand->GetDrcIndex()==-1) continue;
70  DoPidMatch(pidcand,prob);
71  }
72 
73 }
74 
76 {
77 // real track were propagated, reconstructed in Barrel DIRC and Parameterized for Fast simulation
78 // Harphool Kumawat h.kumawat@gsi.de
79  Float_t theta = pidcand->GetMomentum().Theta()*180.0/TMath::Pi();
80  Float_t sigma=0.002;
81  if(theta<45.){
82  sigma=0.001*(3.21659E-01 + 7.48416E-02*theta - 9.87561E-05*theta*theta - 2.47129E-06*theta*theta*theta);
83  }else if(theta>45. && theta<90.){
84  sigma=0.001*(3.21659E-01 + 7.48416E-02*theta - 9.87561E-05*theta*theta - 2.47129E-06*theta*theta*theta);
85  }else if(theta>=90.){
86  sigma=0.001*( 2.34224e+01-3.52791e-01*theta+1.64046e-03*theta*theta-7.57365e-07*theta*theta*theta)+0.0003;
87  }
88  sigma=sigma+0.00873*exp(-pidcand->GetMomentum().Mag()/0.4614);
90  // electron
91  {
92  Float_t mass = 0.0005;
93 // Float_t sigma = 0.002;//0.006; changed by Maria Patsyuk (m.patsyuk@gsi.de) on 3. July 2012
94  prob->SetElectronPdf(GetPdf(pidcand->GetDrcThetaC(),pidcand->GetMomentum().Mag(),mass, sigma));
95  }
96 
97  // muon
98  {
99  Float_t mass = 0.106;
100 // Float_t sigma = 0.002;//0.006;
101  prob->SetMuonPdf(GetPdf(pidcand->GetDrcThetaC(),pidcand->GetMomentum().Mag(),mass, sigma));
102  }
103 
104  // pion
105  {
106  Float_t mass = 0.140;
107 // Float_t sigma = 0.002;//0.006;
108  prob->SetPionPdf(GetPdf(pidcand->GetDrcThetaC(),pidcand->GetMomentum().Mag(),mass, sigma));
109  }
110 
111  // kaon
112  {
113  Float_t mass = 0.494;
114 // Float_t sigma = 0.002;//0.005;
115  prob->SetKaonPdf(GetPdf(pidcand->GetDrcThetaC(),pidcand->GetMomentum().Mag(),mass, sigma));
116  }
117 
118  // proton
119  {
120  Float_t mass = 0.938;
121 // Float_t sigma = 0.002;//0.005;
122  prob->SetProtonPdf(GetPdf(pidcand->GetDrcThetaC(),pidcand->GetMomentum().Mag(),mass, sigma));
123  }
124 }
125 
127 {
128  Double_t beta = mom / TMath::Sqrt(mom*mom + mass*mass);
129  if ( (beta>0.) && ((1./1.47/beta)<1.) )
130  {
131  Double_t center = 1./1.47/beta;
132  TF1 *gausPdf = new TF1("gausPdf","gausn",0,1);
133  gausPdf->SetParameter(0,1);
134  gausPdf->SetParameter(1,center);
135  gausPdf->SetParameter(2,sigma);
136  Double_t val = gausPdf->Eval(TMath::Cos(thetaC));
137  delete gausPdf;
138  return val;
139  }
140  else
141  {
142  // FIXME: Don't write Zeros to that pdf!
143  return 0.;
144  }
145 }
146 
147 //_________________________________________________________________
149  //---
150  FairRootManager::Instance()->
151  Register("PidAlgoDrc"+fTrackBranchNamePidHypo,"Pid", fPidChargedProb, kTRUE);
152 }
153 
154 //_________________________________________________________________
156 }
157 
158 //_________________________________________________________________
160  //---
161 }
162 
163 
int fVerbose
Definition: poormantracks.C:24
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
Int_t i
Definition: run_full.C:25
void SetPionPdf(Double_t val)
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
void SetKaonPdf(Double_t val)
Double_t sigma[nsteps]
Definition: dedx_bands.C:65
Double_t mom
Definition: plot_dirc.C:14
static T Cos(const T &x)
Definition: PndCAMath.h:43
Float_t GetDrcThetaC() const
void SetElectronPdf(Double_t val)
void SetMuonPdf(Double_t val)
h_MC_angle SetTitle("MC truth: opening angle of #pi^{0}")
Double_t
TClonesArray * fPidChargedProb
PndPidCandidate TCA for charged particles.
TString name
Double_t GetPdf(Double_t thetaC, Double_t mom, Double_t mass, Double_t sigma)
void SetProtonPdf(Double_t val)
void DoPidMatch(PndPidCandidate *pidcand, PndPidProbability *prob)
virtual void Exec(Option_t *option)
ClassImp(PndAnaContFact)
Int_t GetDrcIndex() const
TVector3 GetMomentum() const
Double_t Pi
TString fTrackBranchNamePidHypo
PndPidProbability TCA for charged particles.
Double_t thetaC
Definition: plot_dirc.C:16
void SetIndex(Int_t idx)