FairRoot/PandaRoot
PndFsmStt.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmStt.cc,v 1.13 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmStt
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 "PndFsmStt.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  _X0=100.0*11.0/_n;
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  _X0=100.0*11.0/_n;
76  _thtMin=_thtMin*M_PI/180.0;
77  _thtMax=_thtMax*M_PI/180.0;
78  //print(std::cout);
79 }
80 
81 //--------------
82 // Destructor --
83 //--------------
84 
86 {
87 }
88 
89 //--------------
90 // Operations --
91 //--------------
92 
95 {
96  PndFsmResponse *result=new PndFsmResponse();
97 
98  result->setDetector(this);
99  bool wasDetected=detected(t);
100  result->setDetected(wasDetected);
101 
102  if (wasDetected && fabs(t->charge())>1e-8)
103  {
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 dEdx 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 
134  result ->setSttdEdx(measdEdx,sig);
135 
136  if (dEdx_e) result->setLHElectron( gauss(measdEdx,dEdx_e,sig) );
137  if (dEdx_mu) result->setLHMuon( gauss(measdEdx,dEdx_mu,sig) );
138  if (dEdx_pi) result->setLHPion( gauss(measdEdx,dEdx_pi,sig) );
139  if (dEdx_K) result->setLHKaon( gauss(measdEdx,dEdx_K,sig) );
140  if (dEdx_p) result->setLHProton(gauss(measdEdx,dEdx_p,sig) );
141  }
142 
143 
144  return result;
145 }
146 
147 double
148 PndFsmStt::gauss(double x, double x0, double s)
149 {
150  return (1.0/(sqrt(2.0*M_PI)*s))*
151  exp(-(x-x0)*(x-x0)/(2.0*s*s));
152 }
153 
154 double
155 PndFsmStt::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::Stt);
207  } else {
208 
209 
210  double theta = t->p4().Theta();
211  //double p = t->p4().Vect().Mag(); //[R.K. 01/2017] unused variable
212  double p_t = t->p4().Vect().Pt();
213  double charge=t->charge();
214 
215  //only charged particles give signal
216  if (fabs(charge)<0.001) return false;
217 
218  // due to track curvature particle doesn't reach barrel
219  double rho = 3.3356 * p_t / _Bfield;
220  if (_rmin>(2*rho)) return false;
221 
222  //due to helix trajectory particle doesn't hit detector (even with dip angle in tht range)
223  double z=2*rho*asin(_rmin/(2*rho))/tan(theta);
224  double polar=atan2(_rmin,z);
225  if (polar<_thtMin || polar>_thtMax) return false;
226 
227  //finally check for efficiency;
228  return ( _rand->Rndm()<=_efficiency);
229 
230  /* double theta = t->p4().Theta();
231  // double p=t->p4().vect().mag();
232  double p_t=t->p4().Vect().Pt();
233  double charge=t->charge();
234  return ( charge!=0.0 && theta>=_thtMin && theta<=_thtMax && p_t>_pmin && _rand->Rndm()<=_efficiency);*/
235  }
236 }
237 
238 
239 double
241 {
242  TLorentzVector p4=t->p4();
243  double theta=p4.Theta();
244  double mom=p4.Vect().Mag();
245  double beta=p4.Beta();
246 
247  double cont1 = 0.3*sqrt(720.0/(_n+4))*_sigXY*mom*sin(theta)/(_Bfield*_Lpath*_Lpath);
248  double cont2 = 5.23e-2/(_Bfield*beta*sqrt(_Lpath*_X0*sin(theta)));
249  double cont3 = _sigTht/tan(theta);
250 
251  return ( sqrt(cont1*cont1 + cont2*cont2 * cont3*cont3) * mom );
252 }
253 
254 double
255 PndFsmStt::dphi(PndFsmTrack *) const // t //[R.K.03/2017] unused variable(s)
256 {
257  return 0.0007; //to be refined
258 }
259 
260 double
261 PndFsmStt::dtheta(PndFsmTrack *) const // t //[R.K.03/2017] unused variable(s)
262 {
263  return 0.0007; //to be refined
264 }
265 
266 void
268 {
269  o <<"Parameters for detector <"<<detName()<<">"<<endl;
270  o <<" _n = "<<_n<<endl;
271  o <<" _sigXY = "<<_sigXY<<endl;
272  o <<" _Bfield = "<<_Bfield<<endl;
273  o <<" _Lpath = "<<_Lpath<<endl;
274  o <<" _X0 = "<<_X0<<endl;
275  o <<" _sigTht = "<<_sigTht<<endl;
276  o <<" _thtMin = "<<_thtMin<<endl;
277  o <<" _thtMax = "<<_thtMax<<endl;
278  o <<" _radiationLength = "<<_radiationLength<<endl;
279  o <<" _pmin = "<<_pmin<<endl;
280  o <<" _rmin = "<<_rmin<<endl;
281  o <<" _dEdxRes = "<<_dEdxRes << " (rel)"<< endl;
282  o <<" _efficiency = "<<_efficiency<<endl;
283 }
284 
285 void
287 {
289  _n = 11.0;
290  _sigXY = 150.0e-6;
291  _Bfield = 2.0;
292  _Lpath = 0.27;
293  _X0 = 0.0;
294  _sigTht = 6.0e-4;
295  _thtMin = 7.765;
296  _thtMax = 159.44;
297  _radiationLength = 0.0;
298  _pmin = 0.0;
299  _rmin = 0.15;
300  _dEdxRes = 0.2; // 20% dEdx resolution
301  _efficiency =1.0;
302 
303 }
304 
305 bool
306 PndFsmStt::setParameter(std::string &name, double value)
307 {
308  // *****************
309  // include here all parameters which should be settable via tcl
310  // *****************
311 
312  bool knownName=true;
313 
314  if (name == "n")
315  _n=value;
316  else
317  if (name == "sigXY")
318  _sigXY=value;
319  else
320  if (name == "Bfield")
321  _Bfield=value;
322  else
323  if (name == "Lpath")
324  _Lpath=value;
325  else
326  if (name == "X0")
327  _X0=value;
328  else
329  if (name == "sigTht")
330  _sigTht=value;
331  else
332  if (name == "thtMin")
333  _thtMin=value;
334  else
335  if (name == "thtMax")
336  _thtMax=value;
337  else
338  if (name == "radiationLength")
339  _radiationLength=value;
340  else
341  if (name == "pmin")
342  _pmin=value;
343  else
344  if (name == "rmin")
345  _rmin=value;
346  else
347  if (name == "dEdxRes")
348  _dEdxRes=value;
349  else
350  if (name == "efficiency")
351  _efficiency=value;
352  else
353  knownName=false;
354 
355  return knownName;
356 }
357 
Double_t x0
Definition: checkhelixhit.C:70
double dphi(PndFsmTrack *t) const
Definition: PndFsmStt.cxx:255
Double_t p
Definition: anasim.C:58
double _efficiency
Definition: PndFsmAbsDet.h:93
void setdphi(double val)
double gauss(double x, double x0, double s)
Definition: PndFsmStt.cxx:148
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
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void setLHMuon(double val)
double dp(PndFsmTrack *t) const
Definition: PndFsmStt.cxx:240
double _radiationLength
Definition: PndFsmStt.h:93
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
bool hitMapResponse(unsigned int)
Definition: PndFsmTrack.h:86
TLorentzVector s
Definition: Pnd2DStar.C:50
Double_t par[3]
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:110
double _rmin
Definition: PndFsmStt.h:92
virtual ~PndFsmStt()
Definition: PndFsmStt.cxx:85
bool setParameter(std::string &name, double value)
Definition: PndFsmStt.cxx:306
double _X0
Definition: PndFsmStt.h:87
double _Lpath
Definition: PndFsmStt.h:86
double charge()
Definition: PndFsmTrack.h:75
Double_t mom
Definition: plot_dirc.C:14
void parseParameterList(ArgList &par)
double dtheta(PndFsmTrack *t) const
Definition: PndFsmStt.cxx:261
void initParameters()
Definition: PndFsmStt.cxx:286
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
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
TLorentzVector p4()
Definition: PndFsmTrack.h:72
double _sigXY
Definition: PndFsmStt.h:84
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
double compdEdx(double p, double M)
Definition: PndFsmStt.cxx:155
const std::string & detName()
Definition: PndFsmAbsDet.h:74
void setdp(double val)
bool detected(PndFsmTrack *t) const
Definition: PndFsmStt.cxx:203
Double_t z
double _Bfield
Definition: PndFsmStt.h:85
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TString name
bool hitMapValid()
Definition: PndFsmTrack.h:85
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
double X
Definition: anaLmdDigi.C:68
double _thtMax
Definition: PndFsmStt.h:90
Double_t x
double _thtMin
Definition: PndFsmStt.h:89
double _n
Definition: PndFsmStt.h:83
std::string _detName
Definition: PndFsmAbsDet.h:92
double Z
Definition: anaLmdDigi.C:68
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
double _dEdxRes
Definition: PndFsmStt.h:94
double _pmin
Definition: PndFsmStt.h:91
void setdtheta(double val)
void setSttdEdx(double val, double err=0)
double _sigTht
Definition: PndFsmStt.h:88
virtual PndFsmResponse * respond(PndFsmTrack *t)
Definition: PndFsmStt.cxx:94
void setLHPion(double val)
void print(std::ostream &o)
Definition: PndFsmStt.cxx:267