FairRoot/PandaRoot
PndFsmEmcPid.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmEmcBarrel.cc,v 1.11 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmEmcPid
7 //
8 // Implementation of the PID info for EMCs
9 //
10 // This software was developed for the PANDA collaboration. If you
11 // use all or part of it, please give an appropriate acknowledgement.
12 //
13 // Author List:
14 // Klaus Goetzen Original Author
15 //
16 // Copyright Information:
17 // Copyright (C) 2014 GSI
18 //
19 //------------------------------------------------------------------------
20 
21 //-----------------------
22 // This Class's Header --
23 //-----------------------
24 #include "PndFsmEmcPid.h"
25 
26 //-------------
27 // C Headers --
28 //-------------
29 
30 //---------------
31 // C++ Headers --
32 //---------------
33 #include <math.h>
34 #include <iostream>
35 
36 using std::cout;
37 using std::endl;
38 using std::ostream;
39 using std::string;
40 
41 //-------------------------------
42 // Collaborating Class Headers --
43 //-------------------------------
44 
45 #include "ArgList.h"
46 #include "PndFsmResponse.h"
47 #include "PndFsmTrack.h"
48 #include "TFile.h"
49 
50 //-----------------------------------------------------------------------
51 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
52 //-----------------------------------------------------------------------
53 
54 //----------------
55 // Constructors --
56 //----------------
57 
59 {
61 
62  _thtMin=_thtMin*M_PI/180.0;
63  _thtMax=_thtMax*M_PI/180.0;
64  _phiMin=_phiMin*M_PI/180.0;
65  _phiMax=_phiMax*M_PI/180.0;
67 
68  //print(std::cout);
69 }
70 
72 {
74  //set default parameter values and parses a parameter list
75  //i.e. std::list<std::string> of the form
76  //"a=1" "b=2" "c=3"
77  parseParameterList(par);
78 
79  _thtMin=_thtMin*M_PI/180.0;
80  _thtMax=_thtMax*M_PI/180.0;
81  _phiMin=_phiMin*M_PI/180.0;
82  _phiMax=_phiMax*M_PI/180.0;
84 
85  //print(std::cout);
86 }
87 
88 
89 
90 
91 //--------------
92 // Destructor --
93 //--------------
94 
96 {
97 }
98 
99 //--------------
100 // Operations --
101 //--------------
102 
105 {
106  PndFsmResponse *result=new PndFsmResponse();
107 
108  result->setDetector(this);
109  bool wasDetected=detected(t);
110  result->setDetected(wasDetected);
111 
112  if (wasDetected)
113  {
114  int type = t->pdt();
115  double p = t->p4().Vect().Mag();
116  double charge = t->charge();
117 
118  // convert pdg code in type: 0=e+-, 1=mu+-, 2=pi+-, 3=K+-, 4=p, 5=p-bar
119  int idx;
120  if (abs(type)==11) idx=0;
121  else if (abs(type)==13) idx=1;
122  else if (abs(type)==211) idx=2;
123  else if (abs(type)==321) idx=3;
124  else if (type==2212) idx=4;
125  else idx=5;
126 
127  // histogram max p = 10.0
128  double hpmax = _emcPidPdf[0]->GetXaxis()->GetXmax();
129  if (p>hpmax) p=hpmax;
130 
131  // get the bin corresponding to p
132  int currbin = _emcPidPdf[0]->GetXaxis()->FindBin(p);
133 
134  // get the slice containing the pdf emcecal(p) for true particle type
135  TH1D *hpdf = _emcPidPdf[idx]->ProjectionY("_tmppdf",currbin, currbin);
136 
137  // get a random emcecal value from the true distribution
138  // if distribution is empty, choose 0
139  double xsig = 0.;
140  if (hpdf->Integral()>0) xsig = hpdf->GetRandom();
141 
142  // store the value in the response object
143  result->setEmcEcal(xsig);
144 
145  // find the bin corresponding to the value
146  int xsigbin = hpdf->FindBin(xsig);
147 
148  // get the probability values for the different particle types for this emcecal value; sum needed for normalization
149  double P[6], Psum=0.;
150  // find P up e ... pbar and sum
151  for (int k=0;k<6;++k)
152  {
153  P[k] = _emcPidPdf[k]->GetBinContent(currbin,xsigbin);
154  if (k<4) Psum += P[k];
155  }
156 
157  // add P_p or P_pbar depending on charge of particle
158  if (charge>0) Psum+=P[4];
159  else Psum+=P[5];
160 
161  if (Psum<=0.) Psum=1.;
162 
163  result->setLHElectron(P[0]/Psum);
164  result->setLHMuon(P[1]/Psum);
165  result->setLHPion(P[2]/Psum);
166  result->setLHKaon(P[3]/Psum);
167  if (charge>0) result->setLHProton(P[4]/Psum);
168  else result->setLHProton(P[5]/Psum);
169  }
170 
171  return result;
172 }
173 
174 bool
176 {
177  double theta = t->p4().Theta();
178  double phi = t->p4().Phi();
179  double p=t->p4().Vect().Mag();
180  double pt=t->p4().Pt();
181  double charge=t->charge();
182  //int type=abs(t->pdt()); //[R.K. 01/2017] unused variable
183 
184  return ( fabs(charge)>1e-6 && theta>=_thtMin && theta<=_thtMax && phi>=_phiMin && phi<=_phiMax && p>_pmin && pt> _ptmin && _rand->Rndm()<=_efficiency);
185 }
186 
187 
188 void
190 {
191  o <<"Detector <"<<_detName<<">"<<endl;
192  o <<" _pmin = "<<_pmin<<endl;
193  o <<" _ptmin = "<<_ptmin<<endl;
194  o <<" _thtMin = "<<_thtMin<<endl;
195  o <<" _thtMax = "<<_thtMax<<endl;
196  o <<" _phiMin = "<<_phiMin<<endl;
197  o <<" _phiMax = "<<_phiMax<<endl;
198  o <<" _efficiency = "<<_efficiency<<endl;
199  o <<" _parFileName = "<<_parFileName<<endl;
200 }
201 
202 void
204 {
205  _detName = "EmcPid";
206  _pmin = 0.00;
207  _ptmin = 0.00;
208  _thtMin = 22.0;
209  _thtMax = 140.0;
210  _phiMin=-360.;
211  _phiMax=360.;
212  _efficiency=1.0;
213  _parFileName = "$VMCWORKDIR/fastsim/EmcPidPdf.root";
214 }
215 
216 bool
217 PndFsmEmcPid::setParameter(std::string &name, double value)
218 {
219  // *****************
220  // include here all parameters which should be settable via tcl
221  // *****************
222 
223  bool knownName=true;
224 
225  if (name == "pmin")
226  _pmin=value;
227  else
228  if (name == "ptmin")
229  _ptmin=value;
230  else
231  if (name == "thtMin")
232  _thtMin=value;
233  else
234  if (name == "thtMax")
235  _thtMax=value;
236  else
237  if (name == "phiMin")
238  _phiMin=value;
239  else
240  if (name == "phiMax")
241  _phiMax=value;
242  else
243  if (name == "efficiency")
244  _efficiency=value;
245  else
246  knownName=false;
247 
248  return knownName;
249 }
250 
252 {
253 
254  TFile *f=new TFile(_parFileName.c_str());
255 
256  for (int i=0;i<6;i++)
257  {
258  _emcPidPdf[i]=0;
259  }
260 
261  if (f->IsZombie())
262  {
263  cout <<" -W- (PndFsmEmcPid::readParameters) - file "<<_parFileName.c_str()
264  <<" doesn't exist."<<endl;
265  exit(0);
266  }
267  else
268  {
269  _emcPidPdf[0]=(TH2F*)f->Get("hpdf_e");
270  _emcPidPdf[1]=(TH2F*)f->Get("hpdf_mu");
271  _emcPidPdf[2]=(TH2F*)f->Get("hpdf_pi");
272  _emcPidPdf[3]=(TH2F*)f->Get("hpdf_k");
273  _emcPidPdf[4]=(TH2F*)f->Get("hpdf_p");
274  _emcPidPdf[5]=(TH2F*)f->Get("hpdf_pb");
275 
276  for (int i=0;i<6;i++) _emcPidPdf[i]->SetDirectory(0);
277 
278  f->Close();
279  }
280  delete f;
281 
282  return true;
283 }
Double_t p
Definition: anasim.C:58
double _efficiency
Definition: PndFsmAbsDet.h:93
void setLHElectron(double val)
Int_t i
Definition: run_full.C:25
void setLHProton(double val)
std::list< std::string > ArgList
Definition: ArgList.h:7
exit(0)
bool setParameter(std::string &name, double value)
void setLHMuon(double val)
std::string _parFileName
Definition: PndFsmEmcPid.h:89
Double_t par[3]
double charge()
Definition: PndFsmTrack.h:75
void parseParameterList(ArgList &par)
void initParameters()
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
double _phiMax
Definition: PndFsmEmcPid.h:87
double _thtMin
Definition: PndFsmEmcPid.h:84
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
TH2F * _emcPidPdf[6]
Definition: PndFsmEmcPid.h:80
double _pmin
Definition: PndFsmEmcPid.h:82
int idx[MAX]
Definition: autocutx.C:38
TLorentzVector p4()
Definition: PndFsmTrack.h:72
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
double _thtMax
Definition: PndFsmEmcPid.h:85
TFile * f
Definition: bump_analys.C:12
bool readParameters()
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
bool detected(PndFsmTrack *t) const
TString name
void setEmcEcal(double val)
void print(std::ostream &o)
GeV c P
std::string _detName
Definition: PndFsmAbsDet.h:92
double _phiMin
Definition: PndFsmEmcPid.h:86
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
double _ptmin
Definition: PndFsmEmcPid.h:83
virtual PndFsmResponse * respond(PndFsmTrack *t)
void setLHPion(double val)
virtual ~PndFsmEmcPid()