FairRoot/PandaRoot
GFMaterialEffects.cxx
Go to the documentation of this file.
1 /* Copyright 2008-2009, Technische Universitaet Muenchen,
2  Authors: Christian Hoeppner & Sebastian Neubert
3 
4  This file is part of GENFIT.
5 
6  GENFIT is free software: you can redistribute it and/or modify
7  it under the terms of the GNU Lesser General Public License as published
8  by the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  GENFIT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU Lesser General Public License for more details.
15 
16  You should have received a copy of the GNU Lesser General Public License
17  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #include"GFMaterialEffects.h"
21 #include<iostream>
22 #include "stdlib.h"
23 
24 #include"TDatabasePDG.h"
25 #include"TGeoMaterial.h"
26 #include"TGeoManager.h"
27 
28 #include"math.h"
29 #include"assert.h"
30 
31 #define DEBUG 0
32 
33 
35 
36 float MeanExcEnergy_get(int Z);
37 float MeanExcEnergy_get(TGeoMaterial*);
38 
39 
41  fEnergyLossBetheBloch(true), fNoiseBetheBloch(true),
42  fNoiseCoulomb(true),
43  fEnergyLossBrems(true), fNoiseBrems(true),
44  me(0.510998910E-3),
45  fstep(0),
46  fbeta(0),
47  fdedx(0),
48  fgamma(0),
49  fgammaSquare(0),
50  fmatDensity(0),
51  fmatZ(0),
52  fmatA(0),
53  fradiationLength(0),
54  fmEE(0),
55  fpdg(0),
56  fcharge(0),
57  fmass(0) {
58 }
59 
61 }
62 
64  if(finstance == NULL) finstance = new GFMaterialEffects();
65  return finstance;
66 }
67 
69  if(finstance != NULL) {
70  delete finstance;
71  finstance = NULL;
72  }
73 }
74 
75 
76 double GFMaterialEffects::effects(const std::vector<TVector3>& points,
77  const std::vector<double>& pointPaths,
78  const double& mom,
79  const int& pdg,
80  const bool& doNoise,
81  TMatrixT<double>* noise,
82  const TMatrixT<double>* jacobian,
83  const TVector3* directionBefore,
84  const TVector3* directionAfter){
85 
86  //assert(points.size()==pointPaths.size());
87  fpdg = pdg;
88 
89  double momLoss=0.;
90 
91  for(unsigned int i=1;i<points.size();++i){
92  TVector3 dir=points.at(i)-points.at(i-1);
93  double dist=dir.Mag();
94  double realPath = pointPaths.at(i);
95 
96  if (dist > 1.E-8) { // do material effects only if distance is not too small
97  dir*=1./dist; //normalize dir
98 
99  double X(0.);
100 
101  gGeoManager->InitTrack(points.at(i-1).X(),points.at(i-1).Y(),points.at(i-1).Z(),
102  dir.X(),dir.Y(),dir.Z());
103 
104  while(X<dist){
105 
106  getParameters();
107 
108  gGeoManager->FindNextBoundaryAndStep(dist-X);
109  fstep = gGeoManager->GetStep();
110 
111 
112 
113  if(fmatZ>1.E-3){ // don't calculate energy loss for vacuum
114 
115  calcBeta(mom);
116 
118  momLoss += realPath/dist * this->energyLossBetheBloch(mom);
119  if (doNoise && fEnergyLossBetheBloch && fNoiseBetheBloch)
120  this->noiseBetheBloch(mom, noise);
121 
122  if (doNoise && fNoiseCoulomb)
123  this->noiseCoulomb(mom, noise, jacobian, directionBefore, directionAfter);
124 
125  if (fEnergyLossBrems)
126  momLoss += realPath/dist * this->energyLossBrems(mom);
127  if (doNoise && fEnergyLossBrems && fNoiseBrems)
128  this->noiseBrems(mom, noise);
129  }
130  X += fstep;
131  }
132  }
133  }
134 
135  return momLoss;
136 }
137 
138 
139 double GFMaterialEffects::stepper(const double& maxDist,
140  const double& posx,
141  const double& posy,
142  const double& posz,
143  const double& dirx,
144  const double& diry,
145  const double& dirz,
146  const double& mom,
147  const int& pdg){
148 
149  static const double maxPloss = .005; // maximum relative momentum loss allowed
150 
151  fpdg = pdg;
152 
153  gGeoManager->InitTrack(posx,posy,posz,dirx,diry,dirz);
154 
155  double X(0.);
156  double dP = 0.;
157  double momLoss = 0.;
158 
159  while(X<maxDist){
160 
161  getParameters();
162 
163  gGeoManager->FindNextBoundaryAndStep(maxDist-X);
164  fstep = gGeoManager->GetStep();
165 
166 
167 
168  if(fmatZ>1.E-3){ // don't calculate energy loss for vacuum
169 
170  calcBeta(mom);
171 
173  momLoss += this->energyLossBetheBloch(mom);
174 
175  if (fEnergyLossBrems)
176  momLoss += this->energyLossBrems(mom);
177  }
178 
179  if(dP + momLoss > mom*maxPloss){
180  double fraction = (mom*maxPloss-dP)/momLoss;
181  dP+=fraction*momLoss;
182  X+=fraction*fstep;
183  break;
184  }
185 
186  dP += momLoss;
187  X += fstep;
188  }
189 
190  return X;
191 }
192 
193 
195  assert(gGeoManager->GetCurrentVolume()->GetMedium()!=NULL);
196  TGeoMaterial * mat = gGeoManager->GetCurrentVolume()->GetMedium()->GetMaterial();
197  fmatDensity = mat->GetDensity();
198  fmatZ = mat->GetZ();
199  fmatA = mat->GetA();
200  fradiationLength = mat->GetRadLen();
201  fmEE = MeanExcEnergy_get(mat);
202 
203  TParticlePDG * part = TDatabasePDG::Instance()->GetParticle(fpdg);
204  fcharge = part->Charge()/(3.);
205  fmass = part->Mass();
206 }
207 
208 
210  fbeta = mom/sqrt(fmass*fmass+mom*mom);
211 
212  //for numerical stability
213  fgammaSquare = 1.-fbeta*fbeta;
214  if(fgammaSquare>1.E-10) fgammaSquare = 1./fgammaSquare;
215  else fgammaSquare = 1.E10;
217 }
218 
219 
220 
221 //---- Energy-loss and Noise calculations -----------------------------------------
222 
224 
225  // calc fdedx, also needed in noiseBetheBloch!
227  double massRatio = me/fmass;
228  double argument = fgammaSquare*fbeta*fbeta*me*1.E3*2./((1.E-6*fmEE) * sqrt(1+2*sqrt(fgammaSquare)*massRatio + massRatio*massRatio));
229  if (argument <= exp(fbeta*fbeta))
230  fdedx = 0.;
231  else{
232  fdedx *= (log(argument)-fbeta*fbeta); // Bethe-Bloch [MeV/cm]
233  fdedx *= 1.E-3; // in GeV/cm, hence 1.e-3
234  if (fdedx<0.) fdedx = 0;
235  }
236 
237  double DE = fstep * fdedx; //always positive
238  double momLoss = sqrt(mom*mom+2.*sqrt(mom*mom+fmass*fmass)*DE+DE*DE) - mom; //always positive
239 
240  //in vacuum it can numerically happen that momLoss becomes a small negative number. A cut-off at 0.01 eV for momentum loss seems reasonable
241  if(fabs(momLoss)<1.E-11)momLoss=1.E-11;
242  return momLoss;
243 }
244 
245 
247  TMatrixT<double>* noise) const{
248 
249 
250  // ENERGY LOSS FLUCTUATIONS; calculate sigma^2(E);
251  double sigma2E = 0.;
252  double zeta = 153.4E3 * fcharge*fcharge/(fbeta*fbeta) * fmatZ/fmatA * fmatDensity * fstep; // eV
253  double Emax = 2.E9*me*fbeta*fbeta*fgammaSquare / (1. + 2.*fgamma*me/fmass + (me/fmass)*(me/fmass) ); // eV
254  double kappa = zeta/Emax;
255 
256  if (kappa > 0.01) { // Vavilov-Gaussian regime
257  sigma2E += zeta*Emax*(1.-fbeta*fbeta/2.); // eV^2
258  }
259  else { // Urban/Landau approximation
260  double alpha = 0.996;
261  double sigmaalpha = 15.76;
262  // calculate number of collisions Nc
263  double I = 16. * pow(fmatZ, 0.9); // eV
264  double f2 = 0.;
265  if (fmatZ > 2.) f2 = 2./fmatZ;
266  double f1 = 1. - f2;
267  double e2 = 10.*fmatZ*fmatZ; // eV
268  double e1 = pow( (I/pow(e2,f2)), 1./f1); // eV
269 
270  double mbbgg2 = 2.E9*fmass*fbeta*fbeta*fgammaSquare; // eV
271  double Sigma1 = fdedx*1.0E9 * f1/e1 * (log(mbbgg2 / e1) - fbeta*fbeta) / (log(mbbgg2 / I) - fbeta*fbeta) * 0.6; // 1/cm
272  double Sigma2 = fdedx*1.0E9 * f2/e2 * (log(mbbgg2 / e2) - fbeta*fbeta) / (log(mbbgg2 / I) - fbeta*fbeta) * 0.6; // 1/cm
273  double Sigma3 = fdedx*1.0E9 * Emax / ( I*(Emax+I)*log((Emax+I)/I) ) * 0.4; // 1/cm
274 
275  double Nc = (Sigma1 + Sigma2 + Sigma3)*fstep;
276 
277  if (Nc > 50.) { // truncated Landau distribution
278  // calculate sigmaalpha (see GEANT3 manual W5013)
279  double RLAMED = -0.422784 - fbeta*fbeta - log(zeta/Emax);
280  double RLAMAX = 0.60715 + 1.1934*RLAMED +(0.67794 + 0.052382*RLAMED)*exp(0.94753+0.74442*RLAMED);
281  // from lambda max to sigmaalpha=sigma (empirical polynomial)
282  if(RLAMAX <= 1010.) {
283  sigmaalpha = 1.975560
284  +9.898841e-02 *RLAMAX
285  -2.828670e-04 *RLAMAX*RLAMAX
286  +5.345406e-07 *pow(RLAMAX,3.)
287  -4.942035e-10 *pow(RLAMAX,4.)
288  +1.729807e-13 *pow(RLAMAX,5.);
289  }
290  else { sigmaalpha = 1.871887E+01 + 1.296254E-02 *RLAMAX; }
291  // alpha=54.6 corresponds to a 0.9996 maximum cut
292  if(sigmaalpha > 54.6) sigmaalpha=54.6;
293  sigma2E += sigmaalpha*sigmaalpha * zeta*zeta; // eV^2
294  }
295  else { // Urban model
296  double Ealpha = I / (1.-(alpha*Emax/(Emax+I))); // eV
297  double meanE32 = I*(Emax+I)/Emax * (Ealpha-I); // eV^2
298  sigma2E += fstep * (Sigma1*e1*e1 + Sigma2*e2*e2 + Sigma3*meanE32); // eV^2
299  }
300  }
301 
302  sigma2E*=1.E-18; // eV -> GeV
303 
304  // update noise matrix
305  (*noise)[6][6] += (mom*mom+fmass*fmass)/pow(mom,6.)*sigma2E;
306 }
307 
308 
310  TMatrixT<double>* noise,
311  const TMatrixT<double>* jacobian,
312  const TVector3* directionBefore,
313  const TVector3* directionAfter) const{
314 
315  // MULTIPLE SCATTERING; calculate sigma^2
316  // PANDA report PV/01-07 eq(43); linear in step length
317  double sigma2 = 225.E-6/(fbeta*fbeta*mom*mom) * fstep/fradiationLength * fmatZ/(fmatZ+1) * log(159.*pow(fmatZ,-1./3.))/log(287.*pow(fmatZ,-0.5)); // sigma^2 = 225E-6/mom^2 * XX0/fbeta^2 * Z/(Z+1) * ln(159*Z^(-1/3))/ln(287*Z^(-1/2)
318 
319  // noiseBefore
320  TMatrixT<double> noiseBefore(7,7);
321 
322  // calculate euler angles theta, psi (so that directionBefore' points in z' direction)
323  double psi = 0;
324  if (fabs((*directionBefore)[1]) < 1E-14) { // numerical case: arctan(+-inf)=+-PI/2
325  if ((*directionBefore)[0] >= 0.) psi = M_PI/2.;
326  else psi = 3.*M_PI/2.;
327  }
328  else {
329  if ((*directionBefore)[1] > 0.) psi = M_PI - atan((*directionBefore)[0]/(*directionBefore)[1]);
330  else psi = -atan((*directionBefore)[0]/(*directionBefore)[1]);
331  }
332  // cache sin and cos
333  double sintheta = sqrt(1-(*directionBefore)[2]*(*directionBefore)[2]); // theta = arccos(directionBefore[2])
334  double costheta = (*directionBefore)[2];
335  double sinpsi = sin(psi);
336  double cospsi = cos(psi);
337 
338  // calculate NoiseBefore Matrix R M R^T
339  double noiseBefore34 = sigma2 * cospsi * sinpsi * sintheta*sintheta; // noiseBefore_ij = noiseBefore_ji
340  double noiseBefore35 = -sigma2 * costheta * sinpsi * sintheta;
341  double noiseBefore45 = sigma2 * costheta * cospsi * sintheta;
342 
343  noiseBefore[3][3] = sigma2 * (cospsi*cospsi + costheta*costheta - costheta*costheta * cospsi*cospsi);
344  noiseBefore[4][3] = noiseBefore34;
345  noiseBefore[5][3] = noiseBefore35;
346 
347  noiseBefore[3][4] = noiseBefore34;
348  noiseBefore[4][4] = sigma2 * (sinpsi*sinpsi + costheta*costheta * cospsi*cospsi);
349  noiseBefore[5][4] = noiseBefore45;
350 
351  noiseBefore[3][5] = noiseBefore35;
352  noiseBefore[4][5] = noiseBefore45;
353  noiseBefore[5][5] = sigma2 * sintheta*sintheta;
354 
355  TMatrixT<double> jacobianT(7,7);
356  jacobianT = (*jacobian);
357  jacobianT.T();
358 
359  noiseBefore = jacobianT*noiseBefore*(*jacobian); //propagate
360 
361  // noiseAfter
362  TMatrixT<double> noiseAfter(7,7);
363 
364  // calculate euler angles theta, psi (so that A' points in z' direction)
365  psi = 0;
366  if (fabs((*directionAfter)[1]) < 1E-14) { // numerical case: arctan(+-inf)=+-PI/2
367  if ((*directionAfter)[0] >= 0.) psi = M_PI/2.;
368  else psi = 3.*M_PI/2.;
369  }
370  else {
371  if ((*directionAfter)[1] > 0.) psi = M_PI - atan((*directionAfter)[0]/(*directionAfter)[1]);
372  else psi = -atan((*directionAfter)[0]/(*directionAfter)[1]);
373  }
374  // cache sin and cos
375  sintheta = sqrt(1-(*directionAfter)[2]*(*directionAfter)[2]); // theta = arccos(directionAfter[2])
376  costheta = (*directionAfter)[2];
377  sinpsi = sin(psi);
378  cospsi = cos(psi);
379 
380  // calculate NoiseAfter Matrix R M R^T
381  double noiseAfter34 = sigma2 * cospsi * sinpsi * sintheta*sintheta; // noiseAfter_ij = noiseAfter_ji
382  double noiseAfter35 = -sigma2 * costheta * sinpsi * sintheta;
383  double noiseAfter45 = sigma2 * costheta * cospsi * sintheta;
384 
385  noiseAfter[3][3] = sigma2 * (cospsi*cospsi + costheta*costheta - costheta*costheta * cospsi*cospsi);
386  noiseAfter[4][3] = noiseAfter34;
387  noiseAfter[5][3] = noiseAfter35;
388 
389  noiseAfter[3][4] = noiseAfter34;
390  noiseAfter[4][4] = sigma2 * (sinpsi*sinpsi + costheta*costheta * cospsi*cospsi);
391  noiseAfter[5][4] = noiseAfter45;
392 
393  noiseAfter[3][5] = noiseAfter35;
394  noiseAfter[4][5] = noiseAfter45;
395  noiseAfter[5][5] = sigma2 * sintheta*sintheta;
396 
397  //calculate mean of noiseBefore and noiseAfter and update noise
398  (*noise) += 0.5*noiseBefore + 0.5*noiseAfter;
399 
400 }
401 
402 
403 double GFMaterialEffects::energyLossBrems(const double& mom) const{
404 
405  if (fabs(fpdg)!=11) return 0; // only for electrons and positrons
406 
407  #if !defined(BETHE)
408  static const double C[101]={ 0.0,-0.960613E-01, 0.631029E-01,-0.142819E-01, 0.150437E-02,-0.733286E-04, 0.131404E-05, 0.859343E-01,-0.529023E-01, 0.131899E-01,-0.159201E-02, 0.926958E-04,-0.208439E-05,-0.684096E+01, 0.370364E+01,-0.786752E+00, 0.822670E-01,-0.424710E-02, 0.867980E-04,-0.200856E+01, 0.129573E+01,-0.306533E+00, 0.343682E-01,-0.185931E-02, 0.392432E-04, 0.127538E+01,-0.515705E+00, 0.820644E-01,-0.641997E-02, 0.245913E-03,-0.365789E-05, 0.115792E+00,-0.463143E-01, 0.725442E-02,-0.556266E-03, 0.208049E-04,-0.300895E-06,-0.271082E-01, 0.173949E-01,-0.452531E-02, 0.569405E-03,-0.344856E-04, 0.803964E-06, 0.419855E-02,-0.277188E-02, 0.737658E-03,-0.939463E-04, 0.569748E-05,-0.131737E-06,-0.318752E-03, 0.215144E-03,-0.579787E-04, 0.737972E-05,-0.441485E-06, 0.994726E-08, 0.938233E-05,-0.651642E-05, 0.177303E-05,-0.224680E-06, 0.132080E-07,-0.288593E-09,-0.245667E-03, 0.833406E-04,-0.129217E-04, 0.915099E-06,-0.247179E-07, 0.147696E-03,-0.498793E-04, 0.402375E-05, 0.989281E-07,-0.133378E-07,-0.737702E-02, 0.333057E-02,-0.553141E-03, 0.402464E-04,-0.107977E-05,-0.641533E-02, 0.290113E-02,-0.477641E-03, 0.342008E-04,-0.900582E-06, 0.574303E-05, 0.908521E-04,-0.256900E-04, 0.239921E-05,-0.741271E-07,-0.341260E-04, 0.971711E-05,-0.172031E-06,-0.119455E-06, 0.704166E-08, 0.341740E-05,-0.775867E-06,-0.653231E-07, 0.225605E-07,-0.114860E-08,-0.119391E-06, 0.194885E-07, 0.588959E-08,-0.127589E-08, 0.608247E-10};
409  static const double xi=2.51, beta=0.99, vl=0.00004;
410  #endif
411  #if defined(BETHE) // no MIGDAL corrections
412  static const double C[101]={ 0.0, 0.834459E-02, 0.443979E-02,-0.101420E-02, 0.963240E-04,-0.409769E-05, 0.642589E-07, 0.464473E-02,-0.290378E-02, 0.547457E-03,-0.426949E-04, 0.137760E-05,-0.131050E-07,-0.547866E-02, 0.156218E-02,-0.167352E-03, 0.101026E-04,-0.427518E-06, 0.949555E-08,-0.406862E-02, 0.208317E-02,-0.374766E-03, 0.317610E-04,-0.130533E-05, 0.211051E-07, 0.158941E-02,-0.385362E-03, 0.315564E-04,-0.734968E-06,-0.230387E-07, 0.971174E-09, 0.467219E-03,-0.154047E-03, 0.202400E-04,-0.132438E-05, 0.431474E-07,-0.559750E-09,-0.220958E-02, 0.100698E-02,-0.596464E-04,-0.124653E-04, 0.142999E-05,-0.394378E-07, 0.477447E-03,-0.184952E-03,-0.152614E-04, 0.848418E-05,-0.736136E-06, 0.190192E-07,-0.552930E-04, 0.209858E-04, 0.290001E-05,-0.133254E-05, 0.116971E-06,-0.309716E-08, 0.212117E-05,-0.103884E-05,-0.110912E-06, 0.655143E-07,-0.613013E-08, 0.169207E-09, 0.301125E-04,-0.461920E-04, 0.871485E-05,-0.622331E-06, 0.151800E-07,-0.478023E-04, 0.247530E-04,-0.381763E-05, 0.232819E-06,-0.494487E-08,-0.336230E-04, 0.223822E-04,-0.384583E-05, 0.252867E-06,-0.572599E-08, 0.105335E-04,-0.567074E-06,-0.216564E-06, 0.237268E-07,-0.658131E-09, 0.282025E-05,-0.671965E-06, 0.565858E-07,-0.193843E-08, 0.211839E-10, 0.157544E-04,-0.304104E-05,-0.624410E-06, 0.120124E-06,-0.457445E-08,-0.188222E-05,-0.407118E-06, 0.375106E-06,-0.466881E-07, 0.158312E-08, 0.945037E-07, 0.564718E-07,-0.319231E-07, 0.371926E-08,-0.123111E-09};
413  static const double xi=2.10, fbeta=1.00, vl=0.001;
414  #endif
415 
416  double BCUT=10000.; // energy up to which soft bremsstrahlung energy loss is calculated
417 
418  double THIGH=100., CHIGH=50.;
419  double dedxBrems=0.;
420 
421  if(BCUT>0.){
422  double T, kc;
423 
424  if(BCUT>=mom) BCUT=mom; // confine BCUT to mom
425 
426  // T=mom, confined to THIGH
427  // kc=BCUT, confined to CHIGH ??
428  if(mom>=THIGH) {
429  T=THIGH;
430  if(BCUT>=THIGH) kc=CHIGH;
431  else kc=BCUT;
432  }
433  else {
434  T=mom;
435  kc=BCUT;
436  }
437 
438  double E=T+me; // total electron energy
439  if(BCUT>T) kc=T;
440 
441  double X=log(T/me);
442  double Y=log(kc/(E*vl));
443 
444  double XX;
445  int K;
446  double S=0., YY=1.;
447 
448  for (unsigned int I=1; I<=2; ++I) {
449  XX=1.;
450  for (unsigned int J=1; J<=6; ++J) {
451  K=6*I+J-6;
452  S=S+C[K]*XX*YY;
453  XX=XX*X;
454  }
455  YY=YY*Y;
456  }
457 
458  for (unsigned int I=3; I<=6; ++I) {
459  XX=1.;
460  for (unsigned int J=1; J<=6; ++J) {
461  K=6*I+J-6;
462  if(Y<=0.) S=S+C[K]*XX*YY;
463  else S=S+C[K+24]*XX*YY;
464  XX=XX*X;
465  }
466  YY=YY*Y;
467  }
468 
469  double SS=0.;
470  YY=1.;
471 
472  for (unsigned int I=1; I<=2; ++I) {
473  XX=1.;
474  for (unsigned int J=1; J<=5; ++J) {
475  K=5*I+J+55;
476  SS=SS+C[K]*XX*YY;
477  XX=XX*X;
478  }
479  YY=YY*Y;
480  }
481 
482  for (unsigned int I=3; I<=5; ++I) {
483  XX=1.;
484  for (unsigned int J=1; J<=5; ++J) {
485  K=5*I+J+55;
486  if(Y<=0.) SS=SS+C[K]*XX*YY;
487  else SS=SS+C[K+15]*XX*YY;
488  XX=XX*X;
489  }
490  YY=YY*Y;
491  }
492 
493  S=S+fmatZ*SS;
494 
495  if(S>0.){
496  double CORR=1.;
497  #if !defined(BETHE)
498  CORR=1./(1.+0.805485E-10*fmatDensity*fmatZ*E*E/(fmatA*kc*kc)); // MIGDAL correction factor
499  #endif
500 
501  double FAC=fmatZ*(fmatZ+xi)*E*E * pow((kc*CORR/T),beta) / (E+me);
502  if(FAC<=0.) return 0.;
503  dedxBrems=FAC*S;
504 
505  double RAT;
506 
507  if(mom>THIGH) {
508  if(BCUT<THIGH) {
509  RAT=BCUT/mom;
510  S=(1.-0.5*RAT+2.*RAT*RAT/9.);
511  RAT=BCUT/T;
512  S=S/(1.-0.5*RAT+2.*RAT*RAT/9.);
513  }
514  else {
515  RAT=BCUT/mom;
516  S=BCUT*(1.-0.5*RAT+2.*RAT*RAT/9.);
517  RAT=kc/T;
518  S=S/(kc*(1.-0.5*RAT+2.*RAT*RAT/9.));
519  }
520  dedxBrems=dedxBrems*S; // GeV barn
521  }
522 
523  dedxBrems = 0.60221367*fmatDensity*dedxBrems/fmatA; // energy loss dE/dx [GeV/cm]
524  }
525  }
526 
527  if (dedxBrems<0.) dedxBrems = 0;
528 
529  double factor=1.; // positron correction factor
530 
531  if (fpdg==-11){
532  static const double AA=7522100., A1=0.415, A3=0.0021, A5=0.00054;
533 
534  double ETA=0.;
535  if(fmatZ>0.) {
536  double X=log(AA*mom/fmatZ*fmatZ);
537  if(X>-8.) {
538  if(X>=+9.) ETA=1.;
539  else {
540  double W=A1*X+A3*pow(X,3.)+A5*pow(X,5.);
541  ETA=0.5+atan(W)/M_PI;
542  }
543  }
544  }
545 
546  double E0;
547 
548  if(ETA<0.0001) factor=1.E-10;
549  else if (ETA>0.9999) factor=1.;
550  else {
551  E0=BCUT/mom;
552  if(E0>1.) E0=1.;
553  if(E0<1.E-8) factor=1.;
554  else factor = ETA * ( 1.-pow(1.-E0, 1./ETA) ) / E0;
555  }
556  }
557 
558  double DE = fstep * factor*dedxBrems; //always positive
559  double momLoss = sqrt(mom*mom+2.*sqrt(mom*mom+fmass*fmass)*DE+DE*DE) - mom; //always positive
560 
561  return momLoss;
562 }
563 
564 
566  TMatrixT<double>* noise) const{
567 
568  if (fabs(fpdg)!=11) return; // only for electrons and positrons
569 
570  double LX = 1.442695*fstep/fradiationLength;
571  double S2B = mom*mom * ( 1./pow(3.,LX) - 1./pow(4.,LX) );
572  double DEDXB = pow(fabs(S2B),0.5);
573  DEDXB = 1.2E9*DEDXB; //eV
574  double sigma2E = DEDXB*DEDXB; //eV^2
575  sigma2E*=1.E-18; // eV -> GeV
576 
577  (*noise)[6][6] += (mom*mom+fmass*fmass)/pow(mom,6.)*sigma2E;
578 }
579 
580 //---------------------------------------------------------------------------------
581 
583 
584 
585 /*
586 Reference for elemental mean excitation energies at:
587 http://physics.nist.gov/PhysRefData/XrayMassCoef/tab1.html
588 */
589 
590 const int MeanExcEnergy_NELEMENTS = 93; // 0 = vacuum, 1 = hydrogen, 92 = uranium
591 const float MeanExcEnergy_vals[] = {1.e15, 19.2, 41.8, 40.0, 63.7, 76.0, 78., 82.0, 95.0, 115.0, 137.0, 149.0, 156.0, 166.0, 173.0, 173.0, 180.0, 174.0, 188.0, 190.0, 191.0, 216.0, 233.0, 245.0, 257.0, 272.0, 286.0, 297.0, 311.0, 322.0, 330.0, 334.0, 350.0, 347.0, 348.0, 343.0, 352.0, 363.0, 366.0, 379.0, 393.0, 417.0, 424.0, 428.0, 441.0, 449.0, 470.0, 470.0, 469.0, 488.0, 488.0, 487.0, 485.0, 491.0, 482.0, 488.0, 491.0, 501.0, 523.0, 535.0, 546.0, 560.0, 574.0, 580.0, 591.0, 614.0, 628.0, 650.0, 658.0, 674.0, 684.0, 694.0, 705.0, 718.0, 727.0, 736.0, 746.0, 757.0, 790.0, 790.0, 800.0, 810.0, 823.0, 823.0, 830.0, 825.0, 794.0, 827.0, 826.0, 841.0, 847.0, 878.0, 890.0};
592 
593 float MeanExcEnergy_get(int Z){
594  assert(Z>=0&&Z<MeanExcEnergy_NELEMENTS);
595  return MeanExcEnergy_vals[Z];
596 }
597 
598 float MeanExcEnergy_get(TGeoMaterial* mat){
599  if(mat->IsMixture()){
600  double logMEE = 0.;
601  double denom = 0.;
602  TGeoMixture *mix = (TGeoMixture*)mat;
603  for(int i=0;i<mix->GetNelements();++i){
604  int index = int(floor((mix->GetZmixt())[i]));
605  logMEE += 1./(mix->GetAmixt())[i]*(mix->GetWmixt())[i]*(mix->GetZmixt())[i]*log(MeanExcEnergy_get(index));
606  denom += (mix->GetWmixt())[i]*(mix->GetZmixt())[i]*1./(mix->GetAmixt())[i];
607  }
608  logMEE/=denom;
609  return exp(logMEE);
610  }
611  else{ // not a mixture
612  int index = int(floor(mat->GetZ()));
613  return MeanExcEnergy_get(index);
614  }
615 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
void noiseCoulomb(const double &mom, TMatrixT< double > *noise, const TMatrixT< double > *jacobian, const TVector3 *directionBefore, const TVector3 *directionAfter) const
calculation of multiple scattering
void noiseBrems(const double &mom, TMatrixT< double > *noise) const
calculation of energy loss straggeling
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
Int_t i
Definition: run_full.C:25
static const double me
Definition: mzparameters.h:12
void calcBeta(double mom)
sets fbeta, fgamma, fgammasquare; must only be used after calling getParameters() ...
TF1 * f1
Definition: reco_analys2.C:50
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Contains stepper and energy loss/noise matrix calculation.
void noiseBetheBloch(const double &mom, TMatrixT< double > *noise) const
calculation of energy loss straggeling
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:110
static void destruct()
ClassImp(GFMaterialEffects) const int MeanExcEnergy_NELEMENTS
static GFMaterialEffects * finstance
Double_t mom
Definition: plot_dirc.C:14
double energyLossBetheBloch(const double &mom)
Returns energy loss.
double Y
Definition: anaLmdDigi.C:68
TGeoManager * gGeoManager
void getParameters()
sets fmatDensity, fmatZ, fmatA, fradiationLength, fmEE, fcharge, fmass;
float MeanExcEnergy_get(int Z)
double energyLossBrems(const double &mom) const
Returns energy loss.
int Pic_FED Eff_lEE C()
TTree * T
Definition: anaLmdReco.C:32
#define W
Definition: createSTT.C:76
PndMixBackgroundEvents * mix
Definition: runrecoMix.C:67
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double stepper(const double &maxDist, const double &posx, const double &posy, const double &posz, const double &dirx, const double &diry, const double &dirz, const double &mom, const int &pdg)
Returns maximum length so that a specified momentum loss will not be exceeded.
double X
Definition: anaLmdDigi.C:68
double E0[3]
TFile * f2
const float MeanExcEnergy_vals[]
double Z
Definition: anaLmdDigi.C:68
double alpha
Definition: f_Init.h:9
static GFMaterialEffects * getInstance()
double noise
double effects(const std::vector< TVector3 > &points, const std::vector< double > &pointPaths, const double &mom, const int &pdg, const bool &doNoise=false, TMatrixT< double > *noise=NULL, const TMatrixT< double > *jacobian=NULL, const TVector3 *directionBefore=NULL, const TVector3 *directionAfter=NULL)
Calculates energy loss in the travelled path, optional calculation of noise matrix.