FairRoot/PandaRoot
PndFsmSttPid.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: PndFsmSttPid.cc,v 1.13 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class PndFsmSttPid
7 //
8 // Implementation of the STT for the FastSim
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) 2006 GSI
18 //
19 //------------------------------------------------------------------------
20 
21 //-----------------------
22 // This Class's Header --
23 //-----------------------
24 #include "PndFsmSttPid.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 
49 //-----------------------------------------------------------------------
50 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
51 //-----------------------------------------------------------------------
52 
53 //----------------
54 // Constructors --
55 //----------------
56 
58 {
60 
61  _thtMin=_thtMin*M_PI/180.0;
62  _thtMax=_thtMax*M_PI/180.0;
63  //print(std::cout);
64 }
65 
67 {
69  //set default parameter values and parses a parameter list
70  //i.e. std::list<std::string> of the form
71  //"a=1" "b=2" "c=3"
72  parseParameterList(par);
73 
74  _thtMin=_thtMin*M_PI/180.0;
75  _thtMax=_thtMax*M_PI/180.0;
76  //print(std::cout);
77 }
78 
79 //--------------
80 // Destructor --
81 //--------------
82 
84 {
85 }
86 
87 //--------------
88 // Operations --
89 //--------------
90 
93 {
94  PndFsmResponse *result=new PndFsmResponse();
95 
96  result->setDetector(this);
97  bool wasDetected=detected(t);
98  result->setDetected(wasDetected);
99 
100  if (wasDetected && fabs(t->charge())>1e-8)
101  {
102  // now the dEdx information
103  TParticlePDG* part = _fdbPDG->GetParticle(t->pdt());
104  double mass = (part) ? part->Mass() : t->p4().M();
105  double p=t->p4().Vect().Mag();
106 
107  // overall resolution for the dEdx measurement
108  double dEdx = compdEdx(p,mass);
109  double sig = _dEdxRes*dEdx;
110 
111  double m_e = _fdbPDG->GetParticle(11)->Mass();
112  double m_mu = _fdbPDG->GetParticle(13)->Mass();
113  double m_pi = _fdbPDG->GetParticle(211)->Mass();
114  double m_K = _fdbPDG->GetParticle(321)->Mass();
115  double m_p = _fdbPDG->GetParticle(2212)->Mass();
116 
117  // compute the expected dEdx values for all particle types
118  // we need these to determine the pdf's (gaussian around nominal dEdx with res sig)
119 
120  double dEdx_e = compdEdx(p,m_e);
121  double dEdx_mu = compdEdx(p,m_mu);
122  double dEdx_pi = compdEdx(p,m_pi);
123  double dEdx_K = compdEdx(p,m_K);
124  double dEdx_p = compdEdx(p,m_p);
125 
126  double measdEdx = _rand->Gaus(dEdx,sig);
127 
128  result ->setSttdEdx(measdEdx,sig);
129 
130  if (dEdx_e) result->setLHElectron( gauss(measdEdx,dEdx_e,sig) );
131  if (dEdx_mu) result->setLHMuon( gauss(measdEdx,dEdx_mu,sig) );
132  if (dEdx_pi) result->setLHPion( gauss(measdEdx,dEdx_pi,sig) );
133  if (dEdx_K) result->setLHKaon( gauss(measdEdx,dEdx_K,sig) );
134  if (dEdx_p) result->setLHProton(gauss(measdEdx,dEdx_p,sig) );
135  }
136 
137 
138  return result;
139 }
140 
141 double
142 PndFsmSttPid::gauss(double x, double x0, double s)
143 {
144  return (1.0/(sqrt(2.0*M_PI)*s))*
145  exp(-(x-x0)*(x-x0)/(2.0*s*s));
146 }
147 
148 double
149 PndFsmSttPid::compdEdx(double p, double M)
150 {
151  double dEdX;
152 
153  p*=1000;
154  M*=1000;
155 
156  const double Z=10;
157  const double A=20;
158  const double z=1;//charge of incident particle in unit of e
159 
160  double beta;
161  beta=p/sqrt(M*M+p*p);//CalculateBeta(KE,M);
162 
163  double gamma;
164  gamma=1./sqrt(1-beta*beta);
165 
166  const double I=10e-6*Z;//0.000188;//MeV
167  const double me=0.511;//Mev/c2
168 
169  double Wmax;
170  Wmax=(2*me*beta*beta*gamma*gamma) / (1 + 2*gamma*me/M + (me/M)*(me/M));
171 
172  //const double C1=0.1535;//MeV cm2/g
173 
174  double X,X0,X1;
175  double kappa=0.307075;
176  X0=0.201;
177  X1=3;
178  X=log10(beta*gamma);
179  double delta;
180  double C,a;
181  C=-5.217;
182  a=0.196;
183 
184  if(X<=X0)
185  delta=0;
186  else if(X<=X1)
187  delta=2*log(10.)*X+C+a*(X1-X)*(X1-X)*(X1-X);
188  else
189  delta=2*log(10.)*X+C;
190 
191  dEdX= ( kappa * (Z/A) * z*z /(beta*beta)) * (log(2*me*beta*beta*gamma*gamma*Wmax / (I*I)) - 2*beta*beta - delta); //-Dshell??
192 
193  return dEdX;
194 }
195 
196 bool
198 {
199  double theta = t->p4().Theta();
200  double p_t=t->p4().Vect().Perp(TVector3(0.,0.,1.));
201  double charge=t->charge();
202 
203  return ( charge!=0.0 && theta>=_thtMin && theta<=_thtMax && p_t>_ptmin && _rand->Rndm()<=_efficiency);
204 }
205 
206 void
208 {
209  o <<"Parameters for detector <"<<detName()<<">"<<endl;
210  o <<" _thtMin = "<<_thtMin<<endl;
211  o <<" _thtMax = "<<_thtMax<<endl;
212  o <<" _ptmin = "<<_ptmin<<endl;
213  o <<" _dEdxRes = "<<_dEdxRes << " (rel)"<< endl;
214  o <<" _efficiency = "<<_efficiency<<endl;
215 }
216 
217 void
219 {
220  _detName = "SttPid";
221  _thtMin = 7.765;
222  _thtMax = 159.44;
223  _ptmin = 0.15;
224  _dEdxRes = 0.2; // 20% dEdx resolution
225  _efficiency = 1.0;
226 
227 }
228 
229 bool
230 PndFsmSttPid::setParameter(std::string &name, double value)
231 {
232  // *****************
233  // include here all parameters which should be settable via tcl
234  // *****************
235 
236  bool knownName=true;
237 
238  if (name == "thtMin")
239  _thtMin=value;
240  else
241  if (name == "thtMax")
242  _thtMax=value;
243  else
244  if (name == "ptmin")
245  _ptmin=value;
246  else
247  if (name == "dEdxRes")
248  _dEdxRes=value;
249  else
250  if (name == "efficiency")
251  _efficiency=value;
252  else
253  knownName=false;
254 
255  return knownName;
256 }
257 
double compdEdx(double p, double M)
Double_t x0
Definition: checkhelixhit.C:70
double _thtMin
Definition: PndFsmSttPid.h:76
Double_t p
Definition: anasim.C:58
double _efficiency
Definition: PndFsmAbsDet.h:93
void setLHElectron(double val)
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
static const double me
Definition: mzparameters.h:12
TDatabasePDG * _fdbPDG
Definition: PndFsmAbsDet.h:95
void setLHProton(double val)
std::list< std::string > ArgList
Definition: ArgList.h:7
bool detected(PndFsmTrack *t) const
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void setLHMuon(double val)
double _thtMax
Definition: PndFsmSttPid.h:77
TLorentzVector s
Definition: Pnd2DStar.C:50
virtual PndFsmResponse * respond(PndFsmTrack *t)
Double_t par[3]
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:110
double _ptmin
Definition: PndFsmSttPid.h:78
double charge()
Definition: PndFsmTrack.h:75
virtual ~PndFsmSttPid()
void parseParameterList(ArgList &par)
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
Double_t dEdX
Definition: anasim.C:58
int Pic_FED Eff_lEE C()
void initParameters()
Int_t a
Definition: anaLmdDigi.C:126
TLorentzVector p4()
Definition: PndFsmTrack.h:72
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
const std::string & detName()
Definition: PndFsmAbsDet.h:74
double _dEdxRes
Definition: PndFsmSttPid.h:79
Double_t z
bool setParameter(std::string &name, double value)
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TString name
double X
Definition: anaLmdDigi.C:68
Double_t x
std::string _detName
Definition: PndFsmAbsDet.h:92
double gauss(double x, double x0, double s)
void print(std::ostream &o)
double Z
Definition: anaLmdDigi.C:68
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
void setSttdEdx(double val, double err=0)
void setLHPion(double val)