Jpp 19.3.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JPHYSICS::JRadiationFunction Class Reference

Fast implementation of class JRadiation. More...

#include <JRadiationFunction.hh>

Inheritance diagram for JPHYSICS::JRadiationFunction:
JPHYSICS::JRadiation

Public Member Functions

 JRadiationFunction (const JRadiation &radiation, const unsigned int number_of_bins, const double Emin, const double Emax)
 Constructor.
 
virtual double TotalCrossSectionEErad (const double E) const override
 Pair production cross section.
 
virtual double TotalCrossSectionGNrad (const double E) const override
 Photo-nuclear cross section.
 
virtual double CalculateACoeff (const double E) const override
 Ionization a parameter.
 
double SigmaEErad (const double E, const double eps) const
 Pair production cross section.
 
double TotalCrossSectionBrems (const double E) const
 Bremsstrahlung cross section.
 
virtual double SigmaGNrad (const double E, const double eps) const
 Photo-nuclear cross section.
 
double EfromBrems (const double E) const
 Bremsstrahlung shower energy.
 
double ThetaRMSfromBrems (const double E, const double v) const
 Get RMS of scattering angle for Bremsstrahlung.
 
double EfromEErad (const double E) const
 Pair production shower energy.
 
double ThetaRMSfromEErad (const double E, const double v) const
 Get RMS of scattering angle for pair production.
 
double EfromGNrad (const double E) const
 Photo-nuclear shower energy.
 
double ThetaRMSfromGNrad (const double E, const double v) const
 Get RMS of scattering angle for photo-nuclear shower.
 

Protected Member Functions

virtual double IntegralofG (const double E, const double eps) const override
 
double GofZEvrho (const double E, const double v, const double r) const
 

Static Protected Member Functions

static double sigmaGammaPparam (const double eps)
 
static double PhiofEofepsparam (const double E, const double eps)
 
static double le ()
 
static double r0 ()
 
static double Astar ()
 
static double B ()
 
static double BP ()
 

Protected Attributes

JFunction1D_t sigmaEE
 
JFunction1D_t sigmaGN
 
JFunction1D_t Acoeff
 
JFunction2D_t integral
 
const double Z
 
const double A
 
const double Dn
 
const double DnP
 
const int steps
 
const double EminBrems
 
const double EminEErad
 
const double EminGNrad
 

Private Types

typedef JTOOLS::JGridPolint1Function1D_t JFunction1D_t
 
typedef JTOOLS::JMultiFunction< JFunction1D_t, JTOOLS::JMapList< JTOOLS::JPolint1FunctionalGridMap > > JFunction2D_t
 

Detailed Description

Fast implementation of class JRadiation.

In this, the methods

are reimplemented using lookup tables.

Definition at line 32 of file JRadiationFunction.hh.

Member Typedef Documentation

◆ JFunction1D_t

◆ JFunction2D_t

Constructor & Destructor Documentation

◆ JRadiationFunction()

JPHYSICS::JRadiationFunction::JRadiationFunction ( const JRadiation & radiation,
const unsigned int number_of_bins,
const double Emin,
const double Emax )
inline

Constructor.

Parameters
radiationJRadiation object
number_of_binsnumber of bins
Eminminimal muon energy [GeV]
Emaxmaximal muon energy [GeV]

Definition at line 48 of file JRadiationFunction.hh.

51 :
52 JRadiation(radiation)
53 {
54 using namespace JTOOLS;
55
56 const double xmin = log(Emin);
57 const double xmax = log(Emax);
58
59
60 integral.configure(make_grid(number_of_bins, xmin, xmax));
61
62 for (JFunction2D_t::iterator i = integral.begin(); i != integral.end(); ++i) {
63
64 const double x = i->getX();
65 const double E = exp(x);
66
67 const double ymin = log(EminEErad/E);
68 const double ymax = 0.0;
69
70 if (ymax > ymin) {
71
72 i->getY().configure(make_grid(number_of_bins, ymin, ymax));
73
74 for (JFunction1D_t::iterator j = i->getY().begin(); j != i->getY().end(); ++j) {
75
76 const double y = j->getX();
77 const double eps = exp(y) * E;
78 const double z = JRadiation::IntegralofG(E, eps);
79
80 j->getY() = z;
81 }
82 }
83 }
84
86 integral.setExceptionHandler(new JFunction1D_t::JDefaultResult(0.0));
87
88
89 sigmaEE.configure(make_grid(number_of_bins, xmin, xmax));
90
91 for (JFunction1D_t::iterator i = sigmaEE.begin(); i != sigmaEE.end(); ++i) {
92
93 const double E = exp(i->getX());
94 const double y = JRadiation::TotalCrossSectionEErad(E);
95
96 i->getY() = y;
97 }
98
99 sigmaEE.compile();
100
101
102 sigmaGN.configure(make_grid(number_of_bins, xmin, xmax));
103
104 for (JFunction1D_t::iterator i = sigmaGN.begin(); i != sigmaGN.end(); ++i) {
105
106 const double E = exp(i->getX());
107 const double y = JRadiation::TotalCrossSectionGNrad(E);
108
109 i->getY() = y;
110 }
111
112 sigmaGN.compile();
113
114
115 Acoeff.configure(make_grid(number_of_bins, xmin, xmax));
116
117 for (JFunction1D_t::iterator i = Acoeff.begin(); i != Acoeff.end(); ++i) {
118
119 const double E = exp(i->getX());
120 const double y = JRadiation::CalculateACoeff(E);
121
122 i->getY() = y;
123 }
124
125 Acoeff.compile();
126 }
virtual double IntegralofG(const double E, const double eps) const
virtual double TotalCrossSectionEErad(const double E) const
Pair production cross section.
virtual double CalculateACoeff(double E) const
Ionization a parameter.
virtual double TotalCrossSectionGNrad(const double E) const
Photo-nuclear cross section.
JRadiation(const double z, const double a, const int integrsteps, const double eminBrems, const double eminEErad, const double eminGNrad)
Constructor.
Definition JRadiation.hh:56
const double EminEErad
void setExceptionHandler(const typename function_type::supervisor_type &supervisor)
Set the supervisor for handling of exceptions.
void compile()
Compilation.
const double xmax
const double xmin
Auxiliary classes and methods for multi-dimensional interpolations and histograms.
int j
Definition JPolint.hh:801
JGrid< JAbscissa_t > make_grid(const int nx, const JAbscissa_t Xmin, const JAbscissa_t Xmax)
Helper method for JGrid.
Definition JGrid.hh:209

Member Function Documentation

◆ TotalCrossSectionEErad()

virtual double JPHYSICS::JRadiationFunction::TotalCrossSectionEErad ( const double E) const
inlineoverridevirtual

Pair production cross section.

Parameters
Emuon energy [GeV]
Returns
cross section [m^2/g]

Reimplemented from JPHYSICS::JRadiation.

Definition at line 135 of file JRadiationFunction.hh.

136 {
137 const double x = log(E);
138
139 if (x >= sigmaEE. begin()->getX() &&
140 x <= sigmaEE.rbegin()->getX()) {
141
142 try {
143 return sigmaEE(x);
144 }
145 catch(std::exception& error) {}
146 }
147
149 }

◆ TotalCrossSectionGNrad()

virtual double JPHYSICS::JRadiationFunction::TotalCrossSectionGNrad ( const double E) const
inlineoverridevirtual

Photo-nuclear cross section.

Parameters
Emuon energy [GeV]
Returns
cross section [m^2/g]

Reimplemented from JPHYSICS::JRadiation.

Definition at line 158 of file JRadiationFunction.hh.

159 {
160 const double x = log(E);
161
162 if (x >= sigmaGN. begin()->getX() &&
163 x <= sigmaGN.rbegin()->getX()) {
164
165 try {
166 return sigmaGN(x);
167 }
168 catch(std::exception& error) {}
169 }
170
172 }

◆ CalculateACoeff()

virtual double JPHYSICS::JRadiationFunction::CalculateACoeff ( const double E) const
inlineoverridevirtual

Ionization a parameter.

Parameters
Emuon energy [GeV]
Returns
ionization coefficient [GeV*m^2/g]

Reimplemented from JPHYSICS::JRadiation.

Definition at line 181 of file JRadiationFunction.hh.

182 {
183 const double x = log(E);
184
185 if (x >= Acoeff. begin()->getX() &&
186 x <= Acoeff.rbegin()->getX()) {
187
188 try {
189 return Acoeff(x);
190 }
191 catch(std::exception& error) {}
192 }
193
195 }

◆ IntegralofG()

virtual double JPHYSICS::JRadiationFunction::IntegralofG ( const double E,
const double eps ) const
inlineoverrideprotectedvirtual

Reimplemented from JPHYSICS::JRadiation.

Definition at line 198 of file JRadiationFunction.hh.

200 {
201 try {
202
203 const double x = log(E);
204 const double y = log(eps/E);
205 const double z = integral(x, y);
206
207 return z;
208 }
209 catch(std::exception& error) {}
210
211 return JRadiation::IntegralofG(E, eps);
212 }

◆ SigmaEErad()

double JPHYSICS::JRadiation::SigmaEErad ( const double E,
const double eps ) const
inlineinherited

Pair production cross section.

Parameters
Emuon energy [GeV]
epsshower energy [GeV]
Returns
cross section [m^2/g]

Definition at line 80 of file JRadiation.hh.

82 {
83 //routine to calculate dsigma/de in cm^2/GeV for radiating a photon of energy eps
84
85 if(eps<4.*MASS_ELECTRON)return 0.;
86 if(eps>E-0.75*exp(1.0)*pow(Z,1./3.)) return 0.;
87 double zeta;
88 if(E<35.*MASS_MUON)
89 {
90 zeta = 0.;
91 }
92 else
93 {
94 zeta = (0.073*log((E/MASS_MUON)/(1.+1.95e-5*pow(Z,2./3.)*E/MASS_MUON))-0.26);
95 if(zeta<0.)
96 {
97 zeta=0.;
98 }
99 else
100 {
101 zeta = zeta/(0.058*log((E/MASS_MUON)/(1.+5.3e-5*pow(Z,1./3.)*E/MASS_MUON))-0.14);
102 }
103 }
104 double integ = IntegralofG(E,eps);
105 //cout<< eps<<" integ "<< integ<<endl;
106 return integ*(4/(3.*M_PI))*(Z*(Z+zeta)/A)*AVOGADRO*pow((ALPHA_ELECTRO_MAGNETIC*r0()),2.)*(1-eps/E)/eps;
107 }
static double r0()
T pow(const T &x, const double y)
Power .
Definition JMath.hh:97
static const double ALPHA_ELECTRO_MAGNETIC
Electro-Magnetic coupling constant.
static const double MASS_ELECTRON
electron mass [GeV]
static const double AVOGADRO
Avogadro's number.

◆ TotalCrossSectionBrems()

double JPHYSICS::JRadiation::TotalCrossSectionBrems ( const double E) const
inlineinherited

Bremsstrahlung cross section.

Parameters
Emuon energy [GeV]
Returns
cross section [m^2/g]

Definition at line 144 of file JRadiation.hh.

145 {
146 // double delta = (MASS_MUON*MASS_MUON*eps)/(2*E*(E-eps));
147 double delta = (MASS_MUON*MASS_MUON)/(2*E);
148 //double v = eps/E;
149 double Phin = log((B()*pow(Z,-1./3.)*(MASS_MUON+delta*(DnP*sqrt(exp(1.0))-2.)))/(DnP*(MASS_ELECTRON+delta*sqrt(exp(1.0))*B()*pow(Z,-1./3.))));
150 if(Phin<0.)Phin=0.;
151 double Phie = log((BP()*pow(Z,-2./3.)*MASS_MUON)/((1+delta*MASS_MUON/(MASS_ELECTRON*MASS_ELECTRON*sqrt(exp(1.0))))*(MASS_ELECTRON+delta*sqrt(exp(1.0))*BP()*pow(Z,-2./3.))));
152 if(Phie<0.)Phie=0.;
153 //if(eps>E/(1+MASS_MUON*MASS_MUON/(2.*MASS_ELECTRON*E)))Phie=0.;
154 double sig = (16./3.)*ALPHA_ELECTRO_MAGNETIC*AVOGADRO*pow((MASS_ELECTRON/MASS_MUON*r0()),2.0)*Z*(Z*Phin+Phie)/(A);
155 double epsint = log((E-MASS_MUON)/EminBrems)-(E-MASS_MUON-EminBrems)/E+0.375*(pow(E-MASS_MUON,2.)-pow(EminBrems,2.0))/pow(E,2.0);
156 return epsint*sig;//"cross section" in m^2/g multiplied by density and inverted gives mean free path
157 }
static double B()
static double BP()
const double EminBrems
static const double MASS_MUON
muon mass [GeV]

◆ SigmaGNrad()

virtual double JPHYSICS::JRadiation::SigmaGNrad ( const double E,
const double eps ) const
inlinevirtualinherited

Photo-nuclear cross section.

Parameters
Emuon energy [GeV]
epsshower energy [GeV]
Returns
cross section [m^2/g]

Definition at line 167 of file JRadiation.hh.

168 {
169
170 //minimum radiated energy is 0.2 GeV, not very critical formulae become invalid for very much lower
171 //result is given in m^2/g GeV
172 if (eps<0.2) return 0.;
173 if (eps>E-MASS_PROTON) return 0.;
174 double Aeff;
175 if (A==1)
176 {Aeff = 1.;}
177 else
178 {Aeff = (0.22*A + 0.78 * pow(A,0.89));}
179 double sigmaGammaPofeps = (49.2 + 11.1 * log(eps) + 151.8/sqrt(eps))*pow(10,-34.);
180 double Psiofeps = ALPHA_ELECTRO_MAGNETIC/PI*Aeff*AVOGADRO/A*sigmaGammaPofeps/eps;
181 double Denom = 1+eps/LAMBDA*(1+LAMBDA/(2*MASS_PROTON)+eps/LAMBDA);
182 double epsoverE = eps/E;
183 double Numerator = E * E * (1 - epsoverE) / (MASS_MUON * MASS_MUON) * (1 + MASS_MUON * MASS_MUON * epsoverE * epsoverE / (LAMBDA * LAMBDA * (1 - epsoverE)));
184 double Factor = 1 - epsoverE + epsoverE * epsoverE / 2 * (1 + 2 * MASS_MUON * MASS_MUON / (LAMBDA * LAMBDA));
185 double PhiofEofeps = epsoverE - 1 + Factor * log (Numerator / Denom);
186 //cout<<"PhiofEofeps "<<PhiofEofeps<<" Psiofeps "<<Psiofeps<<endl;
187 return Psiofeps*PhiofEofeps;
188 }
static const double MASS_PROTON
proton mass [GeV]
Auxiliary data structure to convert (lambda) function to printable object.
Definition JManip.hh:726

◆ EfromBrems()

double JPHYSICS::JRadiation::EfromBrems ( const double E) const
inlineinherited

Bremsstrahlung shower energy.

Parameters
Emuon energy [GeV]
Returns
shower energy [GeV]

Definition at line 222 of file JRadiation.hh.

223 {
224 //double EminBrems(0.01);
225 //generate according to 1/k from minimum energy to E
226 double Er = 0.0;
227 for (int i = 1000; i != 0; --i) {
228 Er = EminBrems*exp(gRandom->Rndm()*log(E/EminBrems));
229 //check on the extra factor (1-v+3/4v^2)
230 if(gRandom->Rndm()<(1.-Er/E+0.75*pow(Er/E,2))) break;
231 }
232 return Er;
233 }

◆ ThetaRMSfromBrems()

double JPHYSICS::JRadiation::ThetaRMSfromBrems ( const double E,
const double v ) const
inlineinherited

Get RMS of scattering angle for Bremsstrahlung.

Parameters
Emuon energy [GeV]
venergy loss fraction
Returns
angle [rad]

Definition at line 243 of file JRadiation.hh.

244 {
245 using namespace std;
246
247 const double precision = 1.0e-3;
248
249 const double k1 = 0.092 * pow(E, -1.0/3.0);
250 const double k2 = 0.052 * pow(E, -1.0) * pow(Z, -1.0/4.0);
251 const double k3 = 0.220 * pow(E, -0.92);
252 const double k4 = 0.260 * pow(E, -0.91);
253
254 double rms = 0.0;
255
256 if (v <= 0.5) {
257
258 rms = max(min(k1*sqrt(v), k2), k3*v);
259
260 } else {
261
262 const double n = 0.81 * sqrt(E) / (sqrt(E) + 1.8);
263
264 double k5 = 0.0;
265
266 for (double vmin = 0.5, vmax = 1.0; ; ) {
267
268 const double v = 0.5 * (vmin + vmax);
269
270 const double y = k4 * pow(v, 1.0 + n) * pow(1.0 - v, -n);
271
272 if (abs(y - 0.2) <= precision) {
273
274 k5 = y * pow(1.0 - v, 1.0/2.0);
275
276 break;
277 }
278
279 if (y < 0.2)
280 vmin = v;
281 else
282 vmax = v;
283 }
284
285 rms = k4 * pow(v, 1.0 + n) * pow(1.0 - v, -n);
286
287 if (rms >= 0.2) {
288 rms = k5 * pow(1.0 - v, -1.0/2.0);
289 }
290 }
291
292 return rms;
293 }
const int n
Definition JPolint.hh:791

◆ EfromEErad()

double JPHYSICS::JRadiation::EfromEErad ( const double E) const
inlineinherited

Pair production shower energy.

Parameters
Emuon energy [GeV]
Returns
shower energy [GeV]

Definition at line 302 of file JRadiation.hh.

303 {
304 //generate according to 1/k from minimum energy to E
305
306 const double eps =0.2;
307 const double IntGmax = IntegralofG(E,eps)*2.0;
308
309 double Er = 0.0;
310 for (int i = 1000; i != 0; --i)
311 {
312 Er = EminEErad*exp(gRandom->Rndm()*log(E/EminEErad));
313 //check on the extra factor, (1-v) and IntofG
314 double factor = (1.-Er/E)*IntegralofG(E,Er)/IntGmax;
315 if(gRandom->Rndm()<factor) break;
316 }
317 return Er;
318 }

◆ ThetaRMSfromEErad()

double JPHYSICS::JRadiation::ThetaRMSfromEErad ( const double E,
const double v ) const
inlineinherited

Get RMS of scattering angle for pair production.

Parameters
Emuon energy [GeV]
venergy loss fraction
Returns
angle [rad]

Definition at line 328 of file JRadiation.hh.

329 {
330 using namespace std;
331
332 const double a = 8.9e-4;
333 const double b = 1.5e-5;
334 const double c = 0.032;
335 const double d = 1.0;
336 const double e = 0.1;
337
338 const double n = -1.0;
339
340 const double u = v - 2.0*MASS_ELECTRON/E;
341
342 if (u > 0.0)
343 return (2.3 + log(E)) * (1.0/E) * pow(1.0 - v, n) * (u*u) * (1.0/(v*v)) *
344 min(a * pow(v, 1.0/4.0) * (1.0 + b*E) + c*v/(v+d), e);
345 else
346 return 0.0;
347 }
const double a

◆ EfromGNrad()

double JPHYSICS::JRadiation::EfromGNrad ( const double E) const
inlineinherited

Photo-nuclear shower energy.

Parameters
Emuon energy [GeV]
Returns
shower energy [GeV]

Definition at line 356 of file JRadiation.hh.

357 {
358 //generate according to 1/k from minimum energy to E
359
360 double cmax = sigmaGammaPparam(EminGNrad);
361 if (cmax < sigmaGammaPparam(E))
362 cmax=sigmaGammaPparam(E);
363
364 double Pmax = PhiofEofepsparam(E,EminGNrad);
365
366 double Er = 0.0;
367 for (int i = 1000; i != 0; --i) {
368 Er = EminGNrad*exp(gRandom->Rndm()*log(E/EminGNrad));
369 const double factor = PhiofEofepsparam(E,Er)*sigmaGammaPparam(Er)/(cmax*Pmax);
370 if (gRandom->Rndm() < factor) return Er;
371 }
372
373 return Er;
374 }
static double sigmaGammaPparam(const double eps)
static double PhiofEofepsparam(const double E, const double eps)
const double EminGNrad

◆ ThetaRMSfromGNrad()

double JPHYSICS::JRadiation::ThetaRMSfromGNrad ( const double E,
const double v ) const
inlineinherited

Get RMS of scattering angle for photo-nuclear shower.

Parameters
Emuon energy [GeV]
venergy loss fraction
Returns
angle [rad]

Definition at line 384 of file JRadiation.hh.

385 {
386 return 0.0;
387 }

◆ GofZEvrho()

double JPHYSICS::JRadiation::GofZEvrho ( const double E,
const double v,
const double r ) const
inlineprotectedinherited

Definition at line 427 of file JRadiation.hh.

430 {
431 const double ksi = MASS_MUON*MASS_MUON*v*v*(1-r*r)/(4*MASS_ELECTRON*MASS_ELECTRON*(1-v));//
432 const double b = v*v/(2*(1-v));//
433 const double Be = ((2+r*r)*(1+b)+ksi*(3+r*r))*log(1+1/ksi)+(1-r*r-b)/(1+ksi)-(3+r*r);
434 const double Bm = ((1+r*r)*(1+3*b/2)-(1+2*b)*(1-r*r)/ksi)*log(1+ksi)+ksi*(1-r*r-b)/(1+ksi)+(1+2*b)*(1-r*r);
435 const double Ye = (5-r*r+4*b*(1+r*r))/(2*(1+3*b)*log(3+1/ksi)-r*r-2*b*(2-r*r));
436 const double Ym = (4+r*r+3*b*(1+r*r))/((1+r*r)*(1.5+2*b)*log(3+ksi)+1-1.5*r*r);
437 const double Le = log((Astar()*pow(Z,-1./3.)*sqrt((1+ksi)*(1+Ye)))/
438 (1.+(2.*MASS_ELECTRON*sqrt(exp(1.))*Astar()*pow(Z,-1./3.)*(1+ksi)*(1+Ye))/(E*v*(1-r*r))));
439 const double Lm = log(((MASS_MUON/MASS_ELECTRON)*Astar()*pow(Z,-1./3.)*sqrt((1.+1./ksi)*(1.+Ym)))/
440 (1.+(2.*MASS_ELECTRON*sqrt(exp(1.))*Astar()*pow(Z,-1./3.)*(1+ksi)*(1+Ym))/(E*v*(1-r*r))));
441 double Phie = Be*Le; if(Phie<0.)Phie=0.;
442 double Phim = Bm*Lm; if(Phim<0.)Phim=0.;
443 return Phie+(MASS_ELECTRON*MASS_ELECTRON/(MASS_MUON*MASS_MUON))*Phim;
444 }
static double Astar()

◆ sigmaGammaPparam()

static double JPHYSICS::JRadiation::sigmaGammaPparam ( const double eps)
inlinestaticprotectedinherited

Definition at line 469 of file JRadiation.hh.

470 {
471 return (49.2 + 11.1 * log(eps) + 151.8/sqrt(eps));
472 }

◆ PhiofEofepsparam()

static double JPHYSICS::JRadiation::PhiofEofepsparam ( const double E,
const double eps )
inlinestaticprotectedinherited

Definition at line 474 of file JRadiation.hh.

475 {
476 const double Denom = 1+eps/LAMBDA*(1+LAMBDA/(2*MASS_PROTON)+eps/LAMBDA);
477 const double epsoverE = eps/E;
478 const double Numerator = E * E * (1 - epsoverE) / (MASS_MUON * MASS_MUON) * (1 + MASS_MUON * MASS_MUON * epsoverE * epsoverE /
479 (LAMBDA * LAMBDA * (1 - epsoverE)));
480 const double Factor = 1 - epsoverE + epsoverE * epsoverE / 2 * (1 + 2 * MASS_MUON * MASS_MUON / (LAMBDA * LAMBDA));
481 return (epsoverE - 1 + Factor * log (Numerator / Denom));
482 }

◆ le()

static double JPHYSICS::JRadiation::le ( )
inlinestaticprotectedinherited

Definition at line 485 of file JRadiation.hh.

485{ return 3.8616E-13;}

◆ r0()

static double JPHYSICS::JRadiation::r0 ( )
inlinestaticprotectedinherited

Definition at line 486 of file JRadiation.hh.

486{ return 2.817940e-15; }

◆ Astar()

static double JPHYSICS::JRadiation::Astar ( )
inlinestaticprotectedinherited

Definition at line 487 of file JRadiation.hh.

487{ return 183.0; }

◆ B()

static double JPHYSICS::JRadiation::B ( )
inlinestaticprotectedinherited

Definition at line 488 of file JRadiation.hh.

488{ return 183.0; }

◆ BP()

static double JPHYSICS::JRadiation::BP ( )
inlinestaticprotectedinherited

Definition at line 489 of file JRadiation.hh.

489{ return 1429.0; }

Member Data Documentation

◆ sigmaEE

JFunction1D_t JPHYSICS::JRadiationFunction::sigmaEE
protected

Definition at line 214 of file JRadiationFunction.hh.

◆ sigmaGN

JFunction1D_t JPHYSICS::JRadiationFunction::sigmaGN
protected

Definition at line 215 of file JRadiationFunction.hh.

◆ Acoeff

JFunction1D_t JPHYSICS::JRadiationFunction::Acoeff
protected

Definition at line 216 of file JRadiationFunction.hh.

◆ integral

JFunction2D_t JPHYSICS::JRadiationFunction::integral
protected

Definition at line 217 of file JRadiationFunction.hh.

◆ Z

const double JPHYSICS::JRadiation::Z
protectedinherited

Definition at line 491 of file JRadiation.hh.

◆ A

const double JPHYSICS::JRadiation::A
protectedinherited

Definition at line 492 of file JRadiation.hh.

◆ Dn

const double JPHYSICS::JRadiation::Dn
protectedinherited

Definition at line 493 of file JRadiation.hh.

◆ DnP

const double JPHYSICS::JRadiation::DnP
protectedinherited

Definition at line 494 of file JRadiation.hh.

◆ steps

const int JPHYSICS::JRadiation::steps
protectedinherited

Definition at line 495 of file JRadiation.hh.

◆ EminBrems

const double JPHYSICS::JRadiation::EminBrems
protectedinherited

Definition at line 496 of file JRadiation.hh.

◆ EminEErad

const double JPHYSICS::JRadiation::EminEErad
protectedinherited

Definition at line 497 of file JRadiation.hh.

◆ EminGNrad

const double JPHYSICS::JRadiation::EminGNrad
protectedinherited

Definition at line 498 of file JRadiation.hh.


The documentation for this class was generated from the following file: