EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MaterialEffects.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file MaterialEffects.cc
1 /* Copyright 2008-2014, 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 "MaterialEffects.h"
21 #include "Exception.h"
22 #include "IO.h"
23 
24 #include <stdexcept>
25 #include <string>
26 #include <stdlib.h>
27 #include <math.h>
28 #include <assert.h>
29 
30 #include <TDatabasePDG.h>
31 #include <TMath.h>
32 
33 #include <TH1D.h>
34 #include <TFile.h>
35 
36 
37 namespace genfit {
38 
39 MaterialEffects* MaterialEffects::instance_ = nullptr;
40 
41 
43  noEffects_(false),
44  energyLossBetheBloch_(true), noiseBetheBloch_(true),
45  noiseCoulomb_(true),
46  energyLossBrems_(true), noiseBrems_(true),
47  ignoreBoundariesBetweenEqualMaterials_(true),
48  me_(0.510998910E-3),
49  stepSize_(0),
50  dEdx_(0),
51  E_(0),
52  matDensity_(0),
53  matZ_(0),
54  matA_(0),
55  radiationLength_(0),
56  mEE_(0),
57  pdg_(0),
58  charge_(0),
59  mass_(0),
60  mscModelCode_(0),
61  materialInterface_(nullptr),
62  debugLvl_(0)
63 {
64 }
65 
67 {
68  if (materialInterface_ != nullptr) delete materialInterface_;
69 }
70 
72 {
73  if (instance_ == nullptr) instance_ = new MaterialEffects();
74  return instance_;
75 }
76 
78 {
79  if (instance_ != nullptr) {
80  delete instance_;
81  instance_ = nullptr;
82  }
83 }
84 
86 {
87  if (materialInterface_ != nullptr) {
88  std::string msg("MaterialEffects::initMaterialInterface(): Already initialized! ");
89  std::runtime_error err(msg);
90  }
91  materialInterface_ = matIfc;
92 }
93 
94 
95 
96 void MaterialEffects::setMscModel(const std::string& modelName)
97 {
98  if (modelName == std::string("GEANE")) {
99  mscModelCode_ = 0;
100  } else if (modelName == std::string("Highland")) {
101  mscModelCode_ = 1;
102  } else {// throw exception
103  std::string errorMsg = std::string("There is no MSC model called \"") + modelName + "\". Maybe it is not implemented or you misspelled the model name";
104  Exception exc(errorMsg, __LINE__, __FILE__);
105  exc.setFatal();
106  errorOut << exc.what();
107  throw exc;
108  }
109 }
110 
111 
112 double MaterialEffects::effects(const std::vector<RKStep>& steps,
113  int materialsFXStart,
114  int materialsFXStop,
115  const double& mom,
116  const int& pdg,
117  M7x7* noise)
118 {
119 
120  if (debugLvl_ > 0) {
121  debugOut << " MaterialEffects::effects \n";
122  }
123 
124  /*debugOut << "noEffects_ " << noEffects_ << "\n";
125  debugOut << "energyLossBetheBloch_ " << energyLossBetheBloch_ << "\n";
126  debugOut << "noiseBetheBloch_ " << noiseBetheBloch_ << "\n";
127  debugOut << "noiseCoulomb_ " << noiseCoulomb_ << "\n";
128  debugOut << "energyLossBrems_ " << energyLossBrems_ << "\n";
129  debugOut << "noiseBrems_ " << noiseBrems_ << "\n";*/
130 
131 
132  if (noEffects_) return 0.;
133 
134  if (materialInterface_ == nullptr) {
135  std::string msg("MaterialEffects hasn't been initialized with a correct AbsMaterialInterface pointer!");
136  std::runtime_error err(msg);
137  throw err;
138  }
139 
140  bool doNoise(noise != nullptr);
141 
142  pdg_ = pdg;
144 
145  double momLoss = 0.;
146 
147  for ( std::vector<RKStep>::const_iterator it = steps.begin() + materialsFXStart; it != steps.begin() + materialsFXStop; ++it) { // loop over steps
148 
149  double realPath = it->matStep_.stepSize_;
150  if (fabs(realPath) < 1.E-8) {
151  // do material effects only if distance is not too small
152  continue;
153  }
154 
155  if (debugLvl_ > 0) {
156  debugOut << " calculate matFX ";
157  if (doNoise)
158  debugOut << "and noise";
159  debugOut << " for ";
160  debugOut << "stepSize = " << it->matStep_.stepSize_ << "\t";
161  it->matStep_.material_.Print();
162  }
163 
164  double stepSign(1.);
165  if (realPath < 0)
166  stepSign = -1.;
167  realPath = fabs(realPath);
168  stepSize_ = realPath;
169 
170 
171  const Material& currentMaterial = it->matStep_.material_;
172  matDensity_ = currentMaterial.density;
173  matZ_ = currentMaterial.Z;
174  matA_ = currentMaterial.A;
175  radiationLength_ = currentMaterial.radiationLength;
176  mEE_ = currentMaterial.mEE;
177 
178 
179  if (matZ_ > 1.E-3) { // don't calculate energy loss for vacuum
180 
181  momLoss += momentumLoss(stepSign, mom - momLoss, false);
182 
183  if (doNoise){
184  // get values for the "effective" energy of the RK step E_
185  double p(0), gammaSquare(0), gamma(0), betaSquare(0);
186  this->getMomGammaBeta(E_, p, gammaSquare, gamma, betaSquare);
187  double pSquare = p*p;
188 
190  this->noiseBetheBloch(*noise, p, betaSquare, gamma, gammaSquare);
191 
192  if (noiseCoulomb_)
193  this->noiseCoulomb(*noise, *((M1x3*) &it->state7_[3]), pSquare, betaSquare);
194 
196  this->noiseBrems(*noise, pSquare, betaSquare);
197  } // end doNoise
198 
199  }
200 
201  } // end loop over steps
202 
203  if (momLoss >= mom) {
204  Exception exc("MaterialEffects::effects ==> momLoss >= momentum, aborting extrapolation!",__LINE__,__FILE__);
205  exc.setFatal();
206  throw exc;
207  }
208 
209  return momLoss;
210 }
211 
212 
214  M1x7& state7,
215  const double& mom, // momentum
216  double& relMomLoss, // relative momloss for the step will be added
217  const int& pdg,
218  Material& currentMaterial,
219  StepLimits& limits,
220  bool varField)
221 {
222 
223  static const double maxRelMomLoss = .01; // maximum relative momentum loss allowed
224  static const double Pmin = 4.E-3; // minimum momentum for propagation [GeV]
225  static const double minStep = 1.E-4; // 1 µm
226 
227  // check momentum
228  if(mom < Pmin){
229  std::ostringstream sstream;
230  sstream << "MaterialEffects::stepper ==> momentum too low: " << mom*1000. << " MeV";
231  Exception exc(sstream.str(),__LINE__,__FILE__);
232  exc.setFatal();
233  throw exc;
234  }
235 
236  // Trivial cases
237  if (noEffects_)
238  return;
239 
240  if (materialInterface_ == nullptr) {
241  std::string msg("MaterialEffects hasn't been initialized with a correct AbsMaterialInterface pointer!");
242  std::runtime_error err(msg);
243  throw err;
244  }
245 
246  if (relMomLoss > maxRelMomLoss) {
247  limits.setLimit(stp_momLoss, 0);
248  return;
249  }
250 
251 
252  double sMax = limits.getLowestLimitSignedVal(); // signed
253 
254  if (fabs(sMax) < minStep)
255  return;
256 
257 
258 
259  pdg_ = pdg;
261 
262 
263  // make minStep
264  state7[0] += limits.getStepSign() * minStep * state7[3];
265  state7[1] += limits.getStepSign() * minStep * state7[4];
266  state7[2] += limits.getStepSign() * minStep * state7[5];
267 
268  materialInterface_->initTrack(state7[0], state7[1], state7[2],
269  limits.getStepSign() * state7[3], limits.getStepSign() * state7[4], limits.getStepSign() * state7[5]);
270 
271 
272  currentMaterial = materialInterface_->getMaterialParameters();
273  matDensity_ = currentMaterial.density;
274  matZ_ = currentMaterial.Z;
275  matA_ = currentMaterial.A;
276  radiationLength_ = currentMaterial.radiationLength;
277  mEE_ = currentMaterial.mEE;
278 
279  if (debugLvl_ > 0) {
280  debugOut << " currentMaterial "; currentMaterial.Print();
281  }
282 
283  // limit due to momloss
284  double relMomLossPer_cm(0);
285  stepSize_ = 1.; // set stepsize for momLoss calculation
286 
287  if (matZ_ > 1.E-3) { // don't calculate energy loss for vacuum
288  relMomLossPer_cm = this->momentumLoss(limits.getStepSign(), mom, true) / mom;
289  }
290 
291  double maxStepMomLoss = fabs((maxRelMomLoss - fabs(relMomLoss)) / relMomLossPer_cm); // >= 0
292  limits.setLimit(stp_momLoss, maxStepMomLoss);
293 
294  if (debugLvl_ > 0) {
295  debugOut << " momLoss exceeded after a step of " << maxStepMomLoss
296  << "; relMomLoss up to now = " << relMomLoss << "\n";
297  }
298 
299 
300  // now look for boundaries
301  sMax = limits.getLowestLimitSignedVal();
302 
303  stepSize_ = limits.getStepSign() * minStep;
304  M1x3 SA;
305  double boundaryStep(sMax);
306 
307  for (unsigned int i=0; i<100; ++i) {
308  if (debugLvl_ > 0) {
309  debugOut << " find next boundary\n";
310  }
311  double step = materialInterface_->findNextBoundary(rep, state7, boundaryStep, varField);
312 
313  if (debugLvl_ > 0) {
314  if (step == 0) {
315  debugOut << " materialInterface_ returned a step of 0 \n";
316  }
317  }
318 
319  stepSize_ += step;
320  boundaryStep -= step;
321 
322  if (debugLvl_ > 0) {
323  debugOut << " made a step of " << step << "\n";
324  }
325 
327  break;
328 
329  if (fabs(stepSize_) >= fabs(sMax))
330  break;
331 
332  // propagate with found step to boundary
333  rep->RKPropagate(state7, nullptr, SA, step, varField);
334 
335  // make minStep to cross boundary
336  state7[0] += limits.getStepSign() * minStep * state7[3];
337  state7[1] += limits.getStepSign() * minStep * state7[4];
338  state7[2] += limits.getStepSign() * minStep * state7[5];
339 
340  materialInterface_->initTrack(state7[0], state7[1], state7[2],
341  limits.getStepSign() * state7[3], limits.getStepSign() * state7[4], limits.getStepSign() * state7[5]);
342 
344 
345  if (debugLvl_ > 0) {
346  debugOut << " material after step: "; materialAfter.Print();
347  }
348 
349  if (materialAfter != currentMaterial)
350  break;
351  }
352 
354 
355 
356  relMomLoss += relMomLossPer_cm * limits.getLowestLimitVal();
357 }
358 
359 
361 {
362  TParticlePDG* part = TDatabasePDG::Instance()->GetParticle(pdg_);
363  charge_ = int(part->Charge() / 3.); // We only ever use the square
364  mass_ = part->Mass(); // GeV
365 }
366 
367 
369  double& mom, double& gammaSquare, double& gamma, double& betaSquare) const {
370 
371  if (Energy <= mass_) {
372  Exception exc("MaterialEffects::getMomGammaBeta - Energy <= mass",__LINE__,__FILE__);
373  exc.setFatal();
374  throw exc;
375  }
376  gamma = Energy/mass_;
377  gammaSquare = gamma*gamma;
378  betaSquare = 1.-1./gammaSquare;
379  mom = Energy*sqrt(betaSquare);
380 }
381 
382 
383 
384 //---- Energy-loss and Noise calculations -----------------------------------------
385 
386 double MaterialEffects::momentumLoss(double stepSign, double mom, bool linear)
387 {
388  double E0 = hypot(mom, mass_);
389  double step = stepSize_*stepSign; // signed
390 
391 
392  // calc dEdx_, also needed in noiseBetheBloch!
393  // using fourth order Runge Kutta
394  //k1 = f(t0,y0)
395  //k2 = f(t0 + h/2, y0 + h/2 * k1)
396  //k3 = f(t0 + h/2, y0 + h/2 * k2)
397  //k4 = f(t0 + h, y0 + h * k3)
398 
399  // This means in our case:
400  //dEdx1 = dEdx(x0, E0)
401  //dEdx2 = dEdx(x0 + h/2, E1); E1 = E0 + h/2 * dEdx1
402  //dEdx3 = dEdx(x0 + h/2, E2); E2 = E0 + h/2 * dEdx2
403  //dEdx4 = dEdx(x0 + h, E3); E3 = E0 + h * dEdx3
404 
405  double dEdx1 = dEdx(E0); // dEdx(x0,p0)
406 
407  if (linear) {
408  dEdx_ = dEdx1;
409  }
410  else { // RK4
411  double E1 = E0 - dEdx1*step/2.;
412  double dEdx2 = dEdx(E1); // dEdx(x0 + h/2, E0 + h/2 * dEdx1)
413 
414  double E2 = E0 - dEdx2*step/2.;
415  double dEdx3 = dEdx(E2); // dEdx(x0 + h/2, E0 + h/2 * dEdx2)
416 
417  double E3 = E0 - dEdx3*step;
418  double dEdx4 = dEdx(E3); // dEdx(x0 + h, E0 + h * dEdx3)
419 
420  dEdx_ = (dEdx1 + 2.*dEdx2 + 2.*dEdx3 + dEdx4)/6.;
421  }
422 
423  E_ = E0 - dEdx_*step*0.5;
424 
425  double dE = step*dEdx_; // positive for positive stepSign
426 
427  double momLoss(0);
428 
429  if (E0 - dE <= mass_) {
430  // Step would stop particle (E_kin <= 0).
431  return momLoss = mom;
432  }
433  else momLoss = mom - sqrt(pow(E0 - dE, 2) - mass_*mass_); // momLoss; positive for positive stepSign
434 
435  if (debugLvl_ > 0) {
436  debugOut << " MaterialEffects::momentumLoss: mom = " << mom << "; E0 = " << E0
437  << "; dEdx = " << dEdx_
438  << "; dE = " << dE << "; mass = " << mass_ << "\n";
439  }
440 
441  //assert(momLoss * stepSign >= 0);
442 
443  return momLoss;
444 }
445 
446 
447 double MaterialEffects::dEdx(double Energy) const {
448 
449  double mom(0), gammaSquare(0), gamma(0), betaSquare(0);
450  this->getMomGammaBeta(Energy, mom, gammaSquare, gamma, betaSquare);
451 
452  double result(0);
453 
455  result += dEdxBetheBloch(betaSquare, gamma, gammaSquare);
456 
457  if (energyLossBrems_)
458  result += dEdxBrems(mom);
459 
460  return result;
461 }
462 
463 
464 double MaterialEffects::dEdxBetheBloch(double betaSquare, double gamma, double gammaSquare) const
465 {
466  static const double betaGammaMin(0.05);
467  if (betaSquare*gammaSquare < betaGammaMin*betaGammaMin) {
468  Exception exc("MaterialEffects::dEdxBetheBloch ==> beta*gamma < 0.05, Bethe-Bloch implementation not valid anymore!",__LINE__,__FILE__);
469  exc.setFatal();
470  throw exc;
471  }
472 
473  // calc dEdx_, also needed in noiseBetheBloch!
474  double result( 0.307075 * matZ_ / matA_ * matDensity_ / betaSquare * charge_ * charge_ );
475  double massRatio( me_ / mass_ );
476  double argument( gammaSquare * betaSquare * me_ * 1.E3 * 2. / ((1.E-6 * mEE_) *
477  sqrt(1. + 2. * gamma * massRatio + massRatio * massRatio)) );
478  result *= log(argument) - betaSquare; // Bethe-Bloch [MeV/cm]
479  result *= 1.E-3; // in GeV/cm, hence 1.e-3
480  if (result < 0.) {
481  result = 0;
482  }
483 
484  return result;
485 }
486 
487 
488 void MaterialEffects::noiseBetheBloch(M7x7& noise, double mom, double betaSquare, double gamma, double gammaSquare) const
489 {
490  // Code ported from GEANT 3 (erland.F)
491 
492  // ENERGY LOSS FLUCTUATIONS; calculate sigma^2(E);
493  double sigma2E ( 0. );
494  double zeta ( 153.4E3 * charge_ * charge_ / betaSquare * matZ_ / matA_ * matDensity_ * fabs(stepSize_) ); // eV
495  double Emax ( 2.E9 * me_ * betaSquare * gammaSquare / (1. + 2.*gamma * me_ / mass_ + (me_ / mass_) * (me_ / mass_)) ); // eV
496  double kappa ( zeta / Emax );
497 
498  if (kappa > 0.01) { // Vavilov-Gaussian regime
499  sigma2E += zeta * Emax * (1. - betaSquare / 2.); // eV^2
500  } else { // Urban/Landau approximation
501  // calculate number of collisions Nc
502  double I = 16. * pow(matZ_, 0.9); // eV
503  double f2 = 0.;
504  if (matZ_ > 2.) f2 = 2. / matZ_;
505  double f1 = 1. - f2;
506  double e2 = 10.*matZ_ * matZ_; // eV
507  double e1 = pow((I / pow(e2, f2)), 1. / f1); // eV
508 
509  double mbbgg2 = 2.E9 * mass_ * betaSquare * gammaSquare; // eV
510  double Sigma1 = dEdx_ * 1.0E9 * f1 / e1 * (log(mbbgg2 / e1) - betaSquare) / (log(mbbgg2 / I) - betaSquare) * 0.6; // 1/cm
511  double Sigma2 = dEdx_ * 1.0E9 * f2 / e2 * (log(mbbgg2 / e2) - betaSquare) / (log(mbbgg2 / I) - betaSquare) * 0.6; // 1/cm
512  double Sigma3 = dEdx_ * 1.0E9 * Emax / (I * (Emax + I) * log((Emax + I) / I)) * 0.4; // 1/cm
513 
514  double Nc = (Sigma1 + Sigma2 + Sigma3) * fabs(stepSize_);
515 
516  if (Nc > 50.) { // truncated Landau distribution
517  double sigmaalpha = 15.76;
518  // calculate sigmaalpha (see GEANT3 manual W5013)
519  double RLAMED = -0.422784 - betaSquare - log(zeta / Emax);
520  double RLAMAX = 0.60715 + 1.1934 * RLAMED + (0.67794 + 0.052382 * RLAMED) * exp(0.94753 + 0.74442 * RLAMED);
521  // from lambda max to sigmaalpha=sigma (empirical polynomial)
522  if (RLAMAX <= 1010.) {
523  sigmaalpha = 1.975560
524  + 9.898841e-02 * RLAMAX
525  - 2.828670e-04 * RLAMAX * RLAMAX
526  + 5.345406e-07 * pow(RLAMAX, 3.)
527  - 4.942035e-10 * pow(RLAMAX, 4.)
528  + 1.729807e-13 * pow(RLAMAX, 5.);
529  } else { sigmaalpha = 1.871887E+01 + 1.296254E-02 * RLAMAX; }
530  // alpha=54.6 corresponds to a 0.9996 maximum cut
531  if (sigmaalpha > 54.6) sigmaalpha = 54.6;
532  sigma2E += sigmaalpha * sigmaalpha * zeta * zeta; // eV^2
533  } else { // Urban model
534  static const double alpha = 0.996;
535  double Ealpha = I / (1. - (alpha * Emax / (Emax + I))); // eV
536  double meanE32 = I * (Emax + I) / Emax * (Ealpha - I); // eV^2
537  sigma2E += fabs(stepSize_) * (Sigma1 * e1 * e1 + Sigma2 * e2 * e2 + Sigma3 * meanE32); // eV^2
538  }
539  }
540 
541  sigma2E *= 1.E-18; // eV -> GeV
542 
543  // update noise matrix, using linear error propagation from E to q/p
544  noise[6 * 7 + 6] += charge_*charge_/betaSquare / pow(mom, 4) * sigma2E;
545 }
546 
547 
549  const M1x3& direction, double momSquare, double betaSquare) const
550 {
551 
552  // MULTIPLE SCATTERING; calculate sigma^2
553  double sigma2 = 0;
554  assert(mscModelCode_ == 0 || mscModelCode_ == 1);
555  const double step = fabs(stepSize_);
556  const double step2 = step * step;
557  if (mscModelCode_ == 0) {// PANDA report PV/01-07 eq(43); linear in step length
558  sigma2 = 225.E-6 * charge_ * charge_ / (betaSquare * momSquare) * step / radiationLength_ * matZ_ / (matZ_ + 1) * log(159.*pow(matZ_, -1. / 3.)) / log(287.*pow(matZ_, -0.5)); // sigma^2 = 225E-6*z^2/mom^2 * XX0/beta_^2 * Z/(Z+1) * ln(159*Z^(-1/3))/ln(287*Z^(-1/2)
559 
560  } else if (mscModelCode_ == 1) { //Highland not linear in step length formula taken from PDG book 2011 edition
561  double stepOverRadLength = step / radiationLength_;
562  double logCor = (1 + 0.038 * log(stepOverRadLength));
563  sigma2 = 0.0136 * 0.0136 * charge_ * charge_ / (betaSquare * momSquare) * stepOverRadLength * logCor * logCor;
564  }
565  //assert(sigma2 >= 0.0);
566  sigma2 = (sigma2 > 0.0 ? sigma2 : 0.0);
567  //XXX debugOut << "MaterialEffects::noiseCoulomb the MSC variance is " << sigma2 << std::endl;
568 
569  M7x7 noiseAfter; // will hold the new MSC noise to cause by the current stepSize_ length
570  std::fill(noiseAfter.begin(), noiseAfter.end(), 0);
571 
572  const M1x3& a = direction; // as an abbreviation
573  // This calculates the MSC angular spread in the 7D global
574  // coordinate system. See PDG 2010, Sec. 27.3 for formulae.
575  noiseAfter[0 * 7 + 0] = sigma2 * step2 / 3.0 * (1 - a[0]*a[0]);
576  noiseAfter[1 * 7 + 0] = -sigma2 * step2 / 3.0 * a[0]*a[1];
577  noiseAfter[2 * 7 + 0] = -sigma2 * step2 / 3.0 * a[0]*a[2];
578  noiseAfter[3 * 7 + 0] = sigma2 * step * 0.5 * (1 - a[0]*a[0]);
579  noiseAfter[4 * 7 + 0] = -sigma2 * step * 0.5 * a[0]*a[1];
580  noiseAfter[5 * 7 + 0] = -sigma2 * step * 0.5 * a[0]*a[1];
581  noiseAfter[0 * 7 + 1] = noiseAfter[1 * 7 + 0];
582  noiseAfter[1 * 7 + 1] = sigma2 * step2 / 3.0 * (1 - a[1]*a[1]);
583  noiseAfter[2 * 7 + 1] = -sigma2 * step2 / 3.0 * a[1]*a[2];
584  noiseAfter[3 * 7 + 1] = noiseAfter[4 * 7 + 0]; // Cov(x,a_y) = Cov(y,a_x)
585  noiseAfter[4 * 7 + 1] = sigma2 * step * 0.5 * (1 - a[1] * a[1]);
586  noiseAfter[5 * 7 + 1] = -sigma2 * step * 0.5 * a[1]*a[2];
587  noiseAfter[0 * 7 + 2] = noiseAfter[2 * 7 + 0];
588  noiseAfter[1 * 7 + 2] = noiseAfter[2 * 7 + 1];
589  noiseAfter[2 * 7 + 2] = sigma2 * step2 / 3.0 * (1 - a[2]*a[2]);
590  noiseAfter[3 * 7 + 2] = noiseAfter[5 * 7 + 0]; // Cov(z,a_x) = Cov(x,a_z)
591  noiseAfter[4 * 7 + 2] = noiseAfter[5 * 7 + 1]; // Cov(y,a_z) = Cov(z,a_y)
592  noiseAfter[5 * 7 + 2] = sigma2 * step * 0.5 * (1 - a[2]*a[2]);
593  noiseAfter[0 * 7 + 3] = noiseAfter[3 * 7 + 0];
594  noiseAfter[1 * 7 + 3] = noiseAfter[3 * 7 + 1];
595  noiseAfter[2 * 7 + 3] = noiseAfter[3 * 7 + 2];
596  noiseAfter[3 * 7 + 3] = sigma2 * (1 - a[0]*a[0]);
597  noiseAfter[4 * 7 + 3] = -sigma2 * a[0]*a[1];
598  noiseAfter[5 * 7 + 3] = -sigma2 * a[0]*a[2];
599  noiseAfter[0 * 7 + 4] = noiseAfter[4 * 7 + 0];
600  noiseAfter[1 * 7 + 4] = noiseAfter[4 * 7 + 1];
601  noiseAfter[2 * 7 + 4] = noiseAfter[4 * 7 + 2];
602  noiseAfter[3 * 7 + 4] = noiseAfter[4 * 7 + 3];
603  noiseAfter[4 * 7 + 4] = sigma2 * (1 - a[1]*a[1]);
604  noiseAfter[5 * 7 + 4] = -sigma2 * a[1]*a[2];
605  noiseAfter[0 * 7 + 5] = noiseAfter[5 * 7 + 0];
606  noiseAfter[1 * 7 + 5] = noiseAfter[5 * 7 + 1];
607  noiseAfter[2 * 7 + 5] = noiseAfter[5 * 7 + 2];
608  noiseAfter[3 * 7 + 5] = noiseAfter[5 * 7 + 3];
609  noiseAfter[4 * 7 + 5] = noiseAfter[5 * 7 + 4];
610  noiseAfter[5 * 7 + 5] = sigma2 * (1 - a[2]*a[2]);
611 // debugOut << "new noise\n";
612 // RKTools::printDim(noiseAfter, 7,7);
613  for (unsigned int i = 0; i < 7 * 7; ++i) {
614  noise[i] += noiseAfter[i];
615  }
616 }
617 
618 
619 double MaterialEffects::dEdxBrems(double mom) const
620 {
621 
622  // Code ported from GEANT 3 (gbrele.F)
623 
624  if (abs(pdg_) != 11) return 0; // only for electrons and positrons
625 
626 #if !defined(BETHE)
627  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};
628  static const double xi = 2.51, beta = 0.99, vl = 0.00004;
629 #endif
630 #if defined(BETHE) // no MIGDAL corrections
631  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};
632  static const double xi = 2.10, beta = 1.00, vl = 0.001;
633 #endif
634 
635  double BCUT = 10000.; // energy up to which soft bremsstrahlung energy loss is calculated
636 
637  static const double THIGH = 100., CHIGH = 50.;
638  double dedxBrems = 0.;
639 
640  if (BCUT > 0.) {
641  double T, kc;
642 
643  if (BCUT > mom) BCUT = mom; // confine BCUT to mom_
644 
645  // T=mom_, confined to THIGH
646  // kc=BCUT, confined to CHIGH ??
647  if (mom > THIGH) {
648  T = THIGH;
649  if (BCUT >= THIGH)
650  kc = CHIGH;
651  else
652  kc = BCUT;
653  } else {
654  T = mom;
655  kc = BCUT;
656  }
657 
658  double E = T + me_; // total electron energy
659  if (BCUT > T)
660  kc = T;
661 
662  double X = log(T / me_);
663  double Y = log(kc / (E * vl));
664 
665  double XX;
666  int K;
667  double S = 0., YY = 1.;
668 
669  for (unsigned int I = 1; I <= 2; ++I) {
670  XX = 1.;
671  for (unsigned int J = 1; J <= 6; ++J) {
672  K = 6 * I + J - 6;
673  S += C[K] * XX * YY;
674  XX *= X;
675  }
676  YY *= Y;
677  }
678 
679  for (unsigned int I = 3; I <= 6; ++I) {
680  XX = 1.;
681  for (unsigned int J = 1; J <= 6; ++J) {
682  K = 6 * I + J - 6;
683  if (Y > 0.)
684  K += 24;
685  S += C[K] * XX * YY;
686  XX *= X;
687  }
688  YY *= Y;
689  }
690 
691  double SS = 0.;
692  YY = 1.;
693 
694  for (unsigned int I = 1; I <= 2; ++I) {
695  XX = 1.;
696  for (unsigned int J = 1; J <= 5; ++J) {
697  K = 5 * I + J + 55;
698  SS += C[K] * XX * YY;
699  XX *= X;
700  }
701  YY *= Y;
702  }
703 
704  for (unsigned int I = 3; I <= 5; ++I) {
705  XX = 1.;
706  for (unsigned int J = 1; J <= 5; ++J) {
707  K = 5 * I + J + 55;
708  if (Y > 0.)
709  K += 15;
710  SS += C[K] * XX * YY;
711  XX *= X;
712  }
713  YY *= Y;
714  }
715 
716  S += matZ_ * SS;
717 
718  if (S > 0.) {
719  double CORR = 1.;
720 #if !defined(BETHE)
721  CORR = 1. / (1. + 0.805485E-10 * matDensity_ * matZ_ * E * E / (matA_ * kc * kc)); // MIGDAL correction factor
722 #endif
723 
724  // We use exp(beta * log(...) here because pow(..., beta) is
725  // REALLY slow and we don't need ultimate numerical precision
726  // for this approximation.
727  double FAC = matZ_ * (matZ_ + xi) * E * E / (E + me_);
728  if (beta == 1.) // That is the #ifdef BETHE case
729  FAC *= kc * CORR / T;
730  else
731  FAC *= exp(beta * log(kc * CORR / T));
732  if (FAC <= 0.)
733  return 0.;
734  dedxBrems = FAC * S;
735 
736 
737  if (mom >= THIGH) {
738  double RAT;
739  if (BCUT < THIGH) {
740  RAT = BCUT / mom;
741  S = (1. - 0.5 * RAT + 2.*RAT * RAT / 9.);
742  RAT = BCUT / T;
743  S /= 1. - 0.5 * RAT + 2.*RAT * RAT / 9.;
744  } else {
745  RAT = BCUT / mom;
746  S = BCUT * (1. - 0.5 * RAT + 2.*RAT * RAT / 9.);
747  RAT = kc / T;
748  S /= kc * (1. - 0.5 * RAT + 2.*RAT * RAT / 9.);
749  }
750  dedxBrems *= S; // GeV barn
751  }
752 
753  dedxBrems *= 0.60221367 * matDensity_ / matA_; // energy loss dE/dx [GeV/cm]
754  }
755  }
756 
757  if (dedxBrems < 0.)
758  dedxBrems = 0;
759 
760  double factor = 1.; // positron correction factor
761 
762  if (pdg_ == -11) {
763  static const double AA = 7522100., A1 = 0.415, A3 = 0.0021, A5 = 0.00054;
764 
765  double ETA = 0.;
766  if (matZ_ > 0.) {
767  double X = log(AA * mom / (matZ_ * matZ_));
768  if (X > -8.) {
769  if (X >= +9.) ETA = 1.;
770  else {
771  double W = A1 * X + A3 * pow(X, 3.) + A5 * pow(X, 5.);
772  ETA = 0.5 + atan(W) / M_PI;
773  }
774  }
775  }
776 
777  if (ETA < 0.0001)
778  factor = 1.E-10;
779  else if (ETA > 0.9999)
780  factor = 1.;
781  else {
782  double E0 = BCUT / mom;
783  if (E0 > 1.)
784  E0 = 1.;
785  if (E0 < 1.E-8)
786  factor = 1.;
787  else
788  factor = ETA * (1. - pow(1. - E0, 1. / ETA)) / E0;
789  }
790  }
791 
792  return factor * dedxBrems; //always positive
793 }
794 
795 
796 void MaterialEffects::noiseBrems(M7x7& noise, double momSquare, double betaSquare) const
797 {
798  // Code ported from GEANT 3 (erland.F) and simplified
799  // E \approx p is assumed.
800  // the factor 1.44 is not in the original Bethe-Heitler model.
801  // It seems to be some empirical correction copied over from some other project.
802 
803  if (abs(pdg_) != 11) return; // only for electrons and positrons
804 
805  double minusXOverLn2 = -1.442695 * fabs(stepSize_) / radiationLength_;
806  double sigma2E = 1.44*(pow(3., minusXOverLn2) - pow(4., minusXOverLn2)) * momSquare;
807  sigma2E = std::max(sigma2E, 0.0); // must be positive
808 
809  // update noise matrix, using linear error propagation from E to q/p
810  noise[6 * 7 + 6] += charge_*charge_/betaSquare / pow(momSquare, 2) * sigma2E;
811 }
812 
813 
814 void MaterialEffects::setDebugLvl(unsigned int lvl) {
815  debugLvl_ = lvl;
816  if (materialInterface_ and debugLvl_ > 1)
818 }
819 
820 
822  pdg_ = pdg;
823  this->getParticleParameters();
824 
825  stepSize_ = 1;
826 
827  materialInterface_->initTrack(0, 0, 0, 1, 1, 1);
828  auto currentMaterial = materialInterface_->getMaterialParameters();
829  matDensity_ = currentMaterial.density;
830  matZ_ = currentMaterial.Z;
831  matA_ = currentMaterial.A;
832  radiationLength_ = currentMaterial.radiationLength;
833  mEE_ = currentMaterial.mEE;
834 
835  double minMom = 0.00001;
836  double maxMom = 10000;
837  int nSteps(10000);
838  double logStepSize = (log10(maxMom) - log10(minMom)) / (nSteps-1);
839 
840  TH1D hdEdxBethe("dEdxBethe", "dEdxBethe; log10(mom)", nSteps, log10(minMom), log10(maxMom));
841  TH1D hdEdxBrems("dEdxBrems", "dEdxBrems; log10(mom)", nSteps, log10(minMom), log10(maxMom));
842 
843  for (int i=0; i<nSteps; ++i) {
844  double mom = pow(10., log10(minMom) + i*logStepSize);
845  double E = hypot(mom, mass_);
846 
847  energyLossBrems_ = false;
848  energyLossBetheBloch_ = true;
849 
850  try {
851  hdEdxBethe.Fill(log10(mom), dEdx(E));
852  }
853  catch (...) {
854 
855  }
856 
857 
858  //debugOut<< "E = " << E << "; dEdx = " << dEdx(E) <<"\n";
859 
860  energyLossBrems_ = true;
861  energyLossBetheBloch_ = false;
862  try {
863  hdEdxBrems.Fill(log10(mom), dEdx(E));
864  }
865  catch (...) {
866 
867  }
868  }
869 
870  energyLossBrems_ = true;
871  energyLossBetheBloch_ = true;
872 
873  std::string Result;//string which will contain the result
874  std::stringstream convert; // stringstream used for the conversion
875  convert << pdg;//add the value of Number to the characters in the stream
876  Result = convert.str();//set Result to the content of the stream
877 
878  TFile outfile("dEdx_" + TString(Result) + ".root", "recreate");
879  outfile.cd();
880  hdEdxBethe.Write();
881  hdEdxBrems.Write();
882  outfile.Close();
883 }
884 
885 } /* End of namespace genfit */
886 
887