FairRoot/PandaRoot
PndPidMvdAssociatorTask.cxx
Go to the documentation of this file.
2 #include "PndPidCandidate.h"
3 #include "PndPidProbability.h"
4 #include "PndPidMvdPar.h"
5 #include "FairRootManager.h"
6 #include "TMath.h"
7 #include "TF1.h"
8 #include "Riostream.h"
9 
10 
11 //___________________________________________________________
13  //
14  FairRootManager *fManager =FairRootManager::Instance();
15  fManager->Write();
16 }
17 
18 //___________________________________________________________
20  //---
21  fPidChargedProb = new TClonesArray("PndPidProbability");
23 }
24 
25 //___________________________________________________________
26 PndPidMvdAssociatorTask::PndPidMvdAssociatorTask(const char *name, const char *title):FairTask(name)
27 {
28  //---
29  fPidChargedProb = new TClonesArray("PndPidProbability");
31  SetTitle(title);
32  }
33 
34 //___________________________________________________________
36 
37  std::cout << "InitStatus PndPidMvdAssociatorTask::Init()" << std::endl;
38 
39  FairRootManager *fManager =FairRootManager::Instance();
40 
41  fPidChargedCand = (TClonesArray *)fManager->GetObject("PidChargedCand"+fTrackBranchNamePidHypo);
42  if ( ! fPidChargedCand) {
43  std::cout << "-I- PndPidMvdAssociatorTask::Init: No PndPidCandidate array PidChargedCand there!" << std::endl;
44  return kERROR;
45  }
46 
47  Register();
48 
49  mvdPara = new PndPidMvdPar();
50 
51  std::cout << "-I- PndPidMvdAssociatorTask::Init: Success!" << std::endl;
52 
53  return kSUCCESS;
54 }
55 
56 //______________________________________________________
58  //--
59 }
60 //______________________________________________________
62  if (fPidChargedProb->GetEntriesFast() != 0) fPidChargedProb->Clear();
63  if(fVerbose>1) std::cout << "-I- Start PndPidMvdAssociatorTask. "<<std::endl;
64 
65  // Get the Candidates
66  for(Int_t i=0; i<fPidChargedCand->GetEntriesFast(); i++)
67  {
69  PndPidProbability* prob = new((*fPidChargedProb)[i]) PndPidProbability(1.,1.,1.,1.,1.,i);// initializes with equal probability
70  if (pidcand->GetMvdDEDX()==0) continue;
71  DoPidMatch(pidcand,prob);
72  }
73 
74 }
75 
77 {
78  Double_t CanMpv, CanSigma;
79 
80  //Electron
81  CanMpv=mvdPara->GetElectronMpv(pidcand->GetMomentum().Mag());
82  CanSigma=mvdPara->GetElectronSigma(pidcand->GetMomentum().Mag());
83  prob->SetElectronPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
84  //cout << "ele:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetElectronPdf() << endl;
85 
86  //Proton moment below 0.6
87  if (pidcand->GetMomentum().Mag() < 0.6) {
88  CanMpv=mvdPara->GetProtonLowMpv(pidcand->GetMomentum().Mag());
89  CanSigma=mvdPara->GetProtonLowSigma(pidcand->GetMomentum().Mag());
90  prob->SetProtonPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
91  //cout << "proton:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetProtonPdf() << endl;
92  }
93  else {
94  //Proton momenta above 0.6
95  CanMpv=mvdPara->GetProtonHighMpv(pidcand->GetMomentum().Mag());
96  CanSigma=mvdPara->GetProtonHighSigma(pidcand->GetMomentum().Mag());
97  prob->SetProtonPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
98  //cout << "proton:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetProtonPdf() << endl;
99  }
100  //Pion
101  CanMpv=mvdPara->GetPionMpv(pidcand->GetMomentum().Mag());
102  CanSigma=mvdPara->GetPionSigma(pidcand->GetMomentum().Mag());
103  prob->SetPionPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
104  //cout << "pion:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetPionPdf() << endl;
105 
106  //Muon
107  CanMpv=mvdPara->GetMuonMpv(pidcand->GetMomentum().Mag());
108  CanSigma=mvdPara->GetMuonSigma(pidcand->GetMomentum().Mag());
109  prob->SetMuonPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
110  //cout << "muon:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetMuonPdf() << endl;
111 
112  //Kaon
113  CanMpv=mvdPara->GetKaonMpv(pidcand->GetMomentum().Mag());
114  CanSigma=mvdPara->GetKaonSigma(pidcand->GetMomentum().Mag());
115  prob->SetKaonPdf(GetPdf(pidcand->GetMvdDEDX(), CanMpv, CanSigma));
116  //cout << "kaon:\t" << pidcand->GetMomentum().Mag() << "\t" << pidcand->GetMvdDEDX() << "\t" << CanMpv << "\t" << CanSigma << "\t" << prob->GetKaonPdf() << endl;
117 
118  //prob->NormalizeTo(1.);
119 }
120 
122 {
123  TF1 *landauPdf = new TF1("landauPdf","landaun",0,1);
124  landauPdf->SetParameter(0,1);
125  landauPdf->SetParameter(1,Mpv);
126  landauPdf->SetParameter(2,Sigma);
127  Double_t val = landauPdf->Eval(dedx);
128  delete landauPdf;
129  return val;
130 }
131 
132 
133 
134 //_________________________________________________________________
136  //---
137  FairRootManager::Instance()->
138  Register("PidAlgoMvd"+fTrackBranchNamePidHypo,"Pid", fPidChargedProb, kTRUE);
139 }
140 
141 //_________________________________________________________________
143 }
144 
145 //_________________________________________________________________
147  //---
148 }
149 
150 
Double_t GetKaonSigma(Double_t momentum)
int fVerbose
Definition: poormantracks.C:24
Double_t GetProtonHighSigma(Double_t momentum)
Double_t GetPionSigma(Double_t momentum)
Double_t GetProtonLowMpv(Double_t momentum)
TClonesArray * fPidChargedProb
PndPidCandidate TCA for charged particles.
Int_t i
Definition: run_full.C:25
Double_t GetProtonLowSigma(Double_t momentum)
void SetPionPdf(Double_t val)
Double_t GetElectronMpv(Double_t momentum)
Double_t GetPdf(Double_t dedx, Double_t Mpv, Double_t Sigma)
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
Double_t GetElectronSigma(Double_t momentum)
Double_t GetPionMpv(Double_t momentum)
Double_t GetMuonSigma(Double_t momentum)
void SetKaonPdf(Double_t val)
void SetElectronPdf(Double_t val)
void SetMuonPdf(Double_t val)
Double_t GetMuonMpv(Double_t momentum)
h_MC_angle SetTitle("MC truth: opening angle of #pi^{0}")
Double_t
Float_t GetMvdDEDX() const
Double_t GetProtonHighMpv(Double_t momentum)
TString name
void SetProtonPdf(Double_t val)
Double_t GetKaonMpv(Double_t momentum)
ClassImp(PndAnaContFact)
void DoPidMatch(PndPidCandidate *pidcand, PndPidProbability *prob)
virtual void Exec(Option_t *option)
TString fTrackBranchNamePidHypo
PndPidProbability TCA for charged particles.
TVector3 GetMomentum() const