88 const std::vector<double>& pointPaths,
92 TMatrixT<Double_t>*
noise,
93 const TMatrixT<Double_t>* jacobian,
94 const TVector3* directionBefore,
95 const TVector3* directionAfter){
102 for(
unsigned int i=1;i<points.size();++i){
106 TVector3
dir=points.at(i)-points.at(i-1);
107 double dist=dir.Mag();
108 double realPath = pointPaths.at(i);
119 gGeoManager->InitTrack(points.at(i-1).X(),points.at(i-1).Y(),points.at(i-1).Z(),
120 dir.X(),dir.Y(),dir.Z());
128 gGeoManager->FindNextBoundaryAndStep(dist-X);
129 fstep = gGeoManager->GetStep();
158 this->
noiseCoulomb(mom, noise, jacobian, directionBefore, directionAfter);
192 static const double maxPloss = .005;
194 gGeoManager->InitTrack(posx,posy,posz,dirx,diry,dirz);
208 gGeoManager->FindNextBoundaryAndStep(maxDist-X);
209 fstep = gGeoManager->GetStep();
238 if(dP + momLoss > mom*maxPloss){
239 double fraction = (mom*maxPloss-dP)/momLoss;
240 if ((fraction <= 0.) || (fraction >= 1.))
241 throw GFException(std::string(__func__) +
": invalid fraction", __LINE__, __FILE__).setFatal();
242 dP+=fraction*momLoss;
255 if (!gGeoManager->GetCurrentVolume()->GetMedium())
256 throw GFException(std::string(__func__) +
": no medium", __LINE__, __FILE__).setFatal();
257 TGeoMaterial *
mat = gGeoManager->GetCurrentVolume()->GetMedium()->GetMaterial();
269 TParticlePDG *
part = TDatabasePDG::Instance()->GetParticle(
fpdg);
271 fmass = part->Mass();
297 if (
fmass==0.0)
return(0.0);
316 double momLoss = sqrt(mom*mom+2.*sqrt(mom*mom+
fmass*
fmass)*DE+DE*DE) - mom;
319 if(fabs(momLoss)<1.
E-11)momLoss=1.E-11;
327 TMatrixT<double>*
noise)
const{
334 double kappa = zeta/Emax;
340 double alpha = 0.996;
341 double sigmaalpha = 15.76;
343 double I = 16. * pow(
fmatZ, 0.9);
348 double e1 = pow( (I/pow(e2,f2)), 1./f1);
353 double Sigma3 =
fdedx*1.0E9 * Emax / ( I*(Emax+I)*log((Emax+I)/I) ) * 0.4;
355 double Nc = (Sigma1 + Sigma2 + Sigma3)*
fstep;
359 double RLAMED = -0.422784 -
fbeta*
fbeta - log(zeta/Emax);
360 double RLAMAX = 0.60715 + 1.1934*RLAMED +(0.67794 + 0.052382*RLAMED)*exp(0.94753+0.74442*RLAMED);
362 if(RLAMAX <= 1010.) {
363 sigmaalpha = 1.975560
364 +9.898841e-02 *RLAMAX
365 -2.828670e-04 *RLAMAX*RLAMAX
366 +5.345406e-07 *pow(RLAMAX,3.)
367 -4.942035e-10 *pow(RLAMAX,4.)
368 +1.729807e-13 *pow(RLAMAX,5.);
370 else { sigmaalpha = 1.871887E+01 + 1.296254E-02 *RLAMAX; }
372 if(sigmaalpha > 54.6) sigmaalpha=54.6;
373 sigma2E += sigmaalpha*sigmaalpha * zeta*zeta;
376 double Ealpha = I / (1.-(alpha*Emax/(Emax+I)));
377 double meanE32 = I*(Emax+I)/Emax * (Ealpha-I);
378 sigma2E +=
fstep * (Sigma1*e1*e1 + Sigma2*e2*e2 + Sigma3*meanE32);
385 (*noise)[6][6] += (mom*mom+
fmass*
fmass)/pow(mom,6.)*sigma2E;
390 TMatrixT<double>*
noise,
391 const TMatrixT<double>* jacobian,
392 const TVector3* directionBefore,
393 const TVector3* directionAfter)
const{
400 TMatrixT<double> noiseBefore(7,7);
404 if (fabs((*directionBefore)[1]) < 1
E-14) {
405 if ((*directionBefore)[0] >= 0.) psi = M_PI/2.;
406 else psi = 3.*M_PI/2.;
409 if ((*directionBefore)[1] > 0.) psi = M_PI - atan((*directionBefore)[0]/(*directionBefore)[1]);
410 else psi = -atan((*directionBefore)[0]/(*directionBefore)[1]);
413 double sintheta = sqrt(1-(*directionBefore)[2]*(*directionBefore)[2]);
414 double costheta = (*directionBefore)[2];
415 double sinpsi = sin(psi);
416 double cospsi = cos(psi);
419 double noiseBefore34 = sigma2 * cospsi * sinpsi * sintheta*sintheta;
420 double noiseBefore35 = -sigma2 * costheta * sinpsi * sintheta;
421 double noiseBefore45 = sigma2 * costheta * cospsi * sintheta;
423 noiseBefore[3][3] = sigma2 * (cospsi*cospsi + costheta*costheta - costheta*costheta * cospsi*cospsi);
424 noiseBefore[4][3] = noiseBefore34;
425 noiseBefore[5][3] = noiseBefore35;
427 noiseBefore[3][4] = noiseBefore34;
428 noiseBefore[4][4] = sigma2 * (sinpsi*sinpsi + costheta*costheta * cospsi*cospsi);
429 noiseBefore[5][4] = noiseBefore45;
431 noiseBefore[3][5] = noiseBefore35;
432 noiseBefore[4][5] = noiseBefore45;
433 noiseBefore[5][5] = sigma2 * sintheta*sintheta;
435 TMatrixT<double> jacobianT(7,7);
436 jacobianT = (*jacobian);
439 noiseBefore = jacobianT*noiseBefore*(*jacobian);
442 TMatrixT<double> noiseAfter(7,7);
446 if (fabs((*directionAfter)[1]) < 1
E-14) {
447 if ((*directionAfter)[0] >= 0.) psi = M_PI/2.;
448 else psi = 3.*M_PI/2.;
451 if ((*directionAfter)[1] > 0.) psi = M_PI - atan((*directionAfter)[0]/(*directionAfter)[1]);
452 else psi = -atan((*directionAfter)[0]/(*directionAfter)[1]);
455 sintheta = sqrt(1-(*directionAfter)[2]*(*directionAfter)[2]);
456 costheta = (*directionAfter)[2];
461 double noiseAfter34 = sigma2 * cospsi * sinpsi * sintheta*sintheta;
462 double noiseAfter35 = -sigma2 * costheta * sinpsi * sintheta;
463 double noiseAfter45 = sigma2 * costheta * cospsi * sintheta;
465 noiseAfter[3][3] = sigma2 * (cospsi*cospsi + costheta*costheta - costheta*costheta * cospsi*cospsi);
466 noiseAfter[4][3] = noiseAfter34;
467 noiseAfter[5][3] = noiseAfter35;
469 noiseAfter[3][4] = noiseAfter34;
470 noiseAfter[4][4] = sigma2 * (sinpsi*sinpsi + costheta*costheta * cospsi*cospsi);
471 noiseAfter[5][4] = noiseAfter45;
473 noiseAfter[3][5] = noiseAfter35;
474 noiseAfter[4][5] = noiseAfter45;
475 noiseAfter[5][5] = sigma2 * sintheta*sintheta;
478 (*noise) += 0.5*noiseBefore + 0.5*noiseAfter;
485 if (fabs(
fpdg)!=11)
return 0;
488 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};
489 static const double xi=2.51,
beta=0.99, vl=0.00004;
491 #if defined(BETHE) // no MIGDAL corrections 492 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};
493 static const double xi=2.10,
fbeta=1.00, vl=0.001;
498 double THIGH=100., CHIGH=50.;
504 if(BCUT>=mom) BCUT=mom;
510 if(BCUT>=THIGH) kc=CHIGH;
522 double Y=log(kc/(E*vl));
528 for (
unsigned int I=1; I<=2; ++I) {
530 for (
unsigned int J=1; J<=6; ++J) {
538 for (
unsigned int I=3; I<=6; ++I) {
540 for (
unsigned int J=1; J<=6; ++J) {
542 if(Y<=0.) S=S+C[
K]*XX*YY;
543 else S=S+C[K+24]*XX*YY;
552 for (
unsigned int I=1; I<=2; ++I) {
554 for (
unsigned int J=1; J<=5; ++J) {
562 for (
unsigned int I=3; I<=5; ++I) {
564 for (
unsigned int J=1; J<=5; ++J) {
566 if(Y<=0.) SS=SS+C[
K]*XX*YY;
567 else SS=SS+C[K+15]*XX*YY;
582 if(FAC<=0.)
return 0.;
590 S=(1.-0.5*RAT+2.*RAT*RAT/9.);
592 S=S/(1.-0.5*RAT+2.*RAT*RAT/9.);
596 S=BCUT*(1.-0.5*RAT+2.*RAT*RAT/9.);
598 S=S/(kc*(1.-0.5*RAT+2.*RAT*RAT/9.));
600 dedxBrems=dedxBrems*S;
607 if (dedxBrems<0.) dedxBrems = 0;
612 static const double AA=7522100., A1=0.415, A3=0.0021, A5=0.00054;
620 double W=A1*X+A3*pow(X,3.)+A5*pow(X,5.);
621 ETA=0.5+atan(W)/M_PI;
628 if(ETA<0.0001) factor=1.E-10;
629 else if (ETA>0.9999) factor=1.;
633 if(E0<1.
E-8) factor=1.;
634 else factor = ETA * ( 1.-pow(1.-E0, 1./ETA) ) / E0;
638 double DE =
fstep * factor*dedxBrems;
639 double momLoss = sqrt(mom*mom+2.*sqrt(mom*mom+
fmass*
fmass)*DE+DE*DE) - mom;
646 TMatrixT<double>*
noise)
const{
648 if (fabs(
fpdg)!=11)
return;
651 double S2B = mom*mom * ( 1./pow(3.,LX) - 1./pow(4.,LX) );
652 double DEDXB = pow(fabs(S2B),0.5);
654 double sigma2E = DEDXB*DEDXB;
657 (*noise)[6][6] += (mom*mom+
fmass*
fmass)/pow(mom,6.)*sigma2E;
667 const float MeanExcEnergy_vals[
MeanExcEnergy_NELEMENTS] = {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};
671 throw GFException(std::string(__func__) +
": unsupported Z", __LINE__, __FILE__).
setFatal();
676 if(mat->IsMixture()){
679 TGeoMixture *mix = (TGeoMixture*)mat;
680 for(
int i=0;i<mix->GetNelements();++i){
681 int index = int(floor((mix->GetZmixt())[i]));
682 logMEE += 1./(mix->GetAmixt())[i]*(mix->GetWmixt())[i]*(mix->GetZmixt())[i]*log(
MeanExcEnergy_get(index));
683 denom += (mix->GetWmixt())[i]*(mix->GetZmixt())[i]*1./(mix->GetAmixt())[i];
689 int index = int(floor(mat->GetZ()));
static GFMaterialEffects * getInstance()
void noiseBrems(const double &mom, TMatrixT< double > *noise) const
calculation of energy loss straggeling
const int MeanExcEnergy_NELEMENTS
void noiseCoulomb(const double &mom, TMatrixT< double > *noise, const TMatrixT< double > *jacobian, const TVector3 *directionBefore, const TVector3 *directionAfter) const
calculation of multiple scattering
virtual ~GFMaterialEffects()
double energyLossBrems(const double &mom) const
Returns energy loss.
void noiseBetheBloch(const double &mom, TMatrixT< double > *noise) const
calculation of energy loss straggling
constexpr double kc
Speed of light in vacuum in LArSoft units [cm/ns].
double MeanExcEnergy_get(int Z)
const float MeanExcEnergy_vals[MeanExcEnergy_NELEMENTS]
static GFMaterialEffects * finstance
void calcBeta(double mom)
sets fbeta, fgamma, fgammasquare; must only be used after calling getParameters() ...
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
double effects(const std::vector< TVector3 > &points, const std::vector< double > &pointPaths, const double &mom, const int &pdg, const bool &doNoise=false, TMatrixT< Double_t > *noise=NULL, const TMatrixT< Double_t > *jacobian=NULL, const TVector3 *directionBefore=NULL, const TVector3 *directionAfter=NULL)
Calculates energy loss in the travelled path, optional calculation of noise matrix.
GFException & setFatal(bool b=true)
set fatal flag. if this is true, the fit stops for this current track repr.
bool fEnergyLossBetheBloch
void getParameters()
contains energy loss classes
double energyLossBetheBloch(const double &mom)
Returns energy loss.
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.