FairRoot/PandaRoot
PndFsmMvd.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmMvd.cc,v 1.10 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmMvd
7 //
8 // Implementation of the MVD 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 "PndFsmMvd.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 //-----------------------------------------------------------------------
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  //print(std::cout);
65 }
66 
68 {
70  //set default parameter values and parses a parameter list
71  //i.e. std::list<std::string> of the form
72  //"a=1" "b=2" "c=3"
73  parseParameterList(par);
74 
75  _thtMin=_thtMin*M_PI/180.0;
76  _thtMax=_thtMax*M_PI/180.0;
77  //print(std::cout);
78 }
79 
80 //--------------
81 // Destructor --
82 //--------------
83 
85 {
86 }
87 
88 //--------------
89 // Operations --
90 //--------------
91 
94 {
95  PndFsmResponse *result=new PndFsmResponse();
96 
97  result->setDetector(this);
98  bool wasDetected=detected(t);
99  result->setDetected(wasDetected);
100 
101  if (wasDetected && fabs(t->charge())>1e-8)
102  {
103  result->setdV(_vtxRes, _vtxRes, _vtxRes);
104  result->setdp(dp(t));
105  result->setdphi(dphi(t));
106  result->setdtheta(dtheta(t));
107 
108  // now the dEdx information
109  TParticlePDG* part = _fdbPDG->GetParticle(t->pdt());
110  double mass = (part) ? part->Mass() : t->p4().M();
111  double p=t->p4().Vect().Mag();
112 
113  // overall resolution for the tht_c measurement
114  double dEdx = compdEdx(p,mass);
115  double sig = _dEdxRes*dEdx;
116 
117  double m_e = _fdbPDG->GetParticle(11)->Mass();
118  double m_mu = _fdbPDG->GetParticle(13)->Mass();
119  double m_pi = _fdbPDG->GetParticle(211)->Mass();
120  double m_K = _fdbPDG->GetParticle(321)->Mass();
121  double m_p = _fdbPDG->GetParticle(2212)->Mass();
122 
123  // compute the expected dEdx values for all particle types
124  // we need these to determine the pdf's (gaussian around nominal dEdx with res sig)
125 
126  double dEdx_e = compdEdx(p,m_e);
127  double dEdx_mu = compdEdx(p,m_mu);
128  double dEdx_pi = compdEdx(p,m_pi);
129  double dEdx_K = compdEdx(p,m_K);
130  double dEdx_p = compdEdx(p,m_p);
131 
132  double measdEdx = _rand->Gaus(dEdx,sig);
133  result ->setMvddEdx(measdEdx,sig);
134 
135  if (dEdx_e) result->setLHElectron( gauss(measdEdx,dEdx_e,sig) );
136  if (dEdx_mu) result->setLHMuon( gauss(measdEdx,dEdx_mu,sig) );
137  if (dEdx_pi) result->setLHPion( gauss(measdEdx,dEdx_pi,sig) );
138  if (dEdx_K) result->setLHKaon( gauss(measdEdx,dEdx_K,sig) );
139  if (dEdx_p) result->setLHProton(gauss(measdEdx,dEdx_p,sig) );
140 
141  }
142 
143  return result;
144 }
145 
146 double
147 PndFsmMvd::gauss(double x, double x0, double s)
148 {
149  return (1.0/(sqrt(2.0*M_PI)*s))*
150  exp(-(x-x0)*(x-x0)/(2.0*s*s));
151 }
152 
153 
154 double
155 PndFsmMvd::compdEdx(double p, double M)
156 {
157  double dEdX;
158 
159  p*=1000;
160  M*=1000;
161 
162  const double Z=10;
163  const double A=20;
164  const double z=1;//charge of incident particle in unit of e
165 
166  double beta;
167  beta=p/sqrt(M*M+p*p);//CalculateBeta(KE,M);
168 
169  double gamma;
170  gamma=1./sqrt(1-beta*beta);
171 
172  const double I=10e-6*Z;//0.000188;//MeV
173  const double me=0.511;//Mev/c2
174 
175  double Wmax;
176  Wmax=(2*me*beta*beta*gamma*gamma) / (1 + 2*gamma*me/M + (me/M)*(me/M));
177 
178  //const double C1=0.1535;//MeV cm2/g
179 
180  double X,X0,X1;
181  double kappa=0.307075;
182  X0=0.201;
183  X1=3;
184  X=log10(beta*gamma);
185  double delta;
186  double C,a;
187  C=-5.217;
188  a=0.196;
189 
190  if(X<=X0)
191  delta=0;
192  else if(X<=X1)
193  delta=2*log(10.)*X+C+a*(X1-X)*(X1-X)*(X1-X);
194  else
195  delta=2*log(10.)*X+C;
196 
197  dEdX= ( kappa * (Z/A) * z*z /(beta*beta)) * (log(2*me*beta*beta*gamma*gamma*Wmax / (I*I)) - 2*beta*beta - delta); //-Dshell??
198 
199  return dEdX;
200 }
201 
202 bool
204 {
205  if (t->hitMapValid()) {
206  return t->hitMapResponse(FsmDetEnum::Mvd);
207  } else {
208  double theta = t->p4().Theta();
209  double p=t->p4().Vect().Mag();
210  double charge=t->charge();
211  return ( charge!=0.0 && theta>=_thtMin && theta<=_thtMax && p>_pmin && _rand->Rndm()<=_efficiency);
212  }
213 }
214 
215 double
217 {
218  double p=t->p4().Vect().Mag();
219  double Dp=_pRes*p; //10% smearing
220  return ( Dp ); //to be refined
221 }
222 
223 double
224 PndFsmMvd::dphi(PndFsmTrack *) const // t //[R.K.03/2017] unused variable(s)
225 {
226  double Dphi=_phiRes*M_PI/180.0;
227  return Dphi; //to be refined
228 }
229 
230 double
231 PndFsmMvd::dtheta(PndFsmTrack *) const // t //[R.K.03/2017] unused variable(s)
232 {
233  double dt=_thetaRes*M_PI/180.0;
234  return dt; //to be refined
235 }
236 
237 void
239 {
240  o <<"Parameters for detector <"<<detName()<<">"<<endl;
241  o <<" _thtMin = "<<_thtMin<<endl;
242  o <<" _thtMax = "<<_thtMax<<endl;
243  o <<" _radiationLength = "<<_radiationLength<<endl;
244  o <<" _pmin = "<<_pmin<<endl;
245  o <<" _vtxRes = "<<_vtxRes*1e5<<"um"<<endl;
246  o <<" _phiRes = "<<_phiRes << " degree" << endl;
247  o <<" _thetaRes = "<<_thetaRes <<" degree" << endl;
248  o <<" _dEdxRes = "<<_dEdxRes << " (rel)"<< endl;
249  o <<" _efficiency = "<<_efficiency<<endl;
250 }
251 
252 void
254 {
256  _thtMin = 5.0;
257  _thtMax = 160.0;
258  _radiationLength = 0.0;
259  _pmin = 0.0;
260  _vtxRes = 100e-4; //vtx res = 100 um
261  _pRes =0.1;//10% smearing
262  _phiRes =0.01; //0.01 degree smearing
263  _thetaRes =0.01; //0.01 degree smearing
264  _dEdxRes =0.25; // 15% dEdx resolution
265  _efficiency =1.0;
266 }
267 
268 bool
269 PndFsmMvd::setParameter(std::string &name, double value)
270 {
271  // *****************
272  // include here all parameters which should be settable via tcl
273  // *****************
274 
275  bool knownName=true;
276 
277  if (name == "thtMin")
278  _thtMin=value;
279  else
280  if (name == "thtMax")
281  _thtMax=value;
282  else
283  if (name == "radiationLength")
284  _radiationLength=value;
285  else
286  if (name == "pmin")
287  _pmin=value;
288  else
289  if (name == "vtxRes")
290  _vtxRes=value;
291  else
292  if (name == "pRes")
293  _pRes=value;
294  else
295  if (name == "phiRes")
296  _phiRes=value;
297  else
298  if (name == "thetaRes")
299  _thetaRes=value;
300  else
301  if (name == "dEdxRes")
302  _dEdxRes=value;
303  else
304  if (name == "efficiency")
305  _efficiency=value;
306  else
307  knownName=false;
308 
309  return knownName;
310 }
311 
Double_t x0
Definition: checkhelixhit.C:70
Double_t p
Definition: anasim.C:58
double _efficiency
Definition: PndFsmAbsDet.h:93
void setdV(TVector3 v)
void setdphi(double val)
bool setParameter(std::string &name, double value)
Definition: PndFsmMvd.cxx:269
void setLHElectron(double val)
virtual ~PndFsmMvd()
Definition: PndFsmMvd.cxx:84
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
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
double dp(PndFsmTrack *t) const
Definition: PndFsmMvd.cxx:216
double compdEdx(double p, double M)
Definition: PndFsmMvd.cxx:155
void setLHMuon(double val)
double _vtxRes
Definition: PndFsmMvd.h:87
bool hitMapResponse(unsigned int)
Definition: PndFsmTrack.h:86
TLorentzVector s
Definition: Pnd2DStar.C:50
double _dEdxRes
Definition: PndFsmMvd.h:91
Double_t par[3]
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:110
double charge()
Definition: PndFsmTrack.h:75
double _thetaRes
Definition: PndFsmMvd.h:90
double _thtMax
Definition: PndFsmMvd.h:84
void parseParameterList(ArgList &par)
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
double _phiRes
Definition: PndFsmMvd.h:89
Double_t dEdX
Definition: anasim.C:58
int Pic_FED Eff_lEE C()
Int_t a
Definition: anaLmdDigi.C:126
static const std::string & name(unsigned int t)
Definition: FsmDetTypes.h:26
double dtheta(PndFsmTrack *t) const
Definition: PndFsmMvd.cxx:231
TLorentzVector p4()
Definition: PndFsmTrack.h:72
double _radiationLength
Definition: PndFsmMvd.h:85
void print(std::ostream &o)
Definition: PndFsmMvd.cxx:238
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
double gauss(double x, double x0, double s)
Definition: PndFsmMvd.cxx:147
const std::string & detName()
Definition: PndFsmAbsDet.h:74
void setdp(double val)
Double_t z
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TString name
bool hitMapValid()
Definition: PndFsmTrack.h:85
double X
Definition: anaLmdDigi.C:68
double _pmin
Definition: PndFsmMvd.h:86
Double_t x
std::string _detName
Definition: PndFsmAbsDet.h:92
double dphi(PndFsmTrack *t) const
Definition: PndFsmMvd.cxx:224
double Z
Definition: anaLmdDigi.C:68
void initParameters()
Definition: PndFsmMvd.cxx:253
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
virtual PndFsmResponse * respond(PndFsmTrack *t)
Definition: PndFsmMvd.cxx:93
void setdtheta(double val)
double _thtMin
Definition: PndFsmMvd.h:83
bool detected(PndFsmTrack *t) const
Definition: PndFsmMvd.cxx:203
double _pRes
Definition: PndFsmMvd.h:88
void setMvddEdx(double val, double err=0)
void setLHPion(double val)