Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
GFlashSamplingShowerParameterisation Class Reference

#include <GFlashSamplingShowerParameterisation.hh>

Inheritance diagram for GFlashSamplingShowerParameterisation:
GVFlashShowerParameterisation

Public Member Functions

 GFlashSamplingShowerParameterisation (G4Material *aMat1, G4Material *aMat2, G4double d1, G4double d2, GFlashSamplingShowerTuning *aPar=0)
 
 ~GFlashSamplingShowerParameterisation ()
 
void ComputeRadialParameters (G4double y, G4double Tau)
 
void GenerateLongitudinalProfile (G4double Energy)
 
void ComputeZAX0EFFetc ()
 
G4double IntegrateEneLongitudinal (G4double LongitudinalStep)
 
G4double IntegrateNspLongitudinal (G4double LongitudinalStep)
 
G4double ComputeTau (G4double LongitudinalPosition)
 
void SetMaterial (G4Material *mat1, G4Material *mat2)
 
G4double GeneratePhi ()
 
G4double GenerateRadius (G4int ispot, G4double Energy, G4double LongitudinalPosition)
 
G4double GenerateExponential (G4double Energy)
 
G4double GetAveR99 ()
 
G4double GetAveR90 ()
 
G4double GetAveTmx ()
 
G4double GetAveT99 ()
 
G4double GetAveT90 ()
 
G4double GetNspot ()
 
G4double GetX0 ()
 
G4double GetEc ()
 
G4double GetRm ()
 
G4double ApplySampling (const G4double DEne, const G4double Energy)
 
- Public Member Functions inherited from GVFlashShowerParameterisation
 GVFlashShowerParameterisation ()
 
virtual ~GVFlashShowerParameterisation ()
 
G4double GeneratePhi ()
 
G4double GetEffZ (const G4Material *material)
 
G4double GetEffA (const G4Material *material)
 
G4double gam (G4double x, G4double a) const
 
void PrintMaterial (const G4Material *mat)
 

Additional Inherited Members

- Protected Attributes inherited from GVFlashShowerParameterisation
GVFlashHomoShowerTuningthePar
 
G4double density
 
G4double A
 
G4double Z
 
G4double X0
 
G4double Ec
 
G4double Rm
 
G4double NSpot
 

Detailed Description

Definition at line 50 of file GFlashSamplingShowerParameterisation.hh.

Constructor & Destructor Documentation

GFlashSamplingShowerParameterisation::GFlashSamplingShowerParameterisation ( G4Material aMat1,
G4Material aMat2,
G4double  d1,
G4double  d2,
GFlashSamplingShowerTuning aPar = 0 
)

Definition at line 48 of file GFlashSamplingShowerParameterisation.cc.

References ComputeZAX0EFFetc(), GFlashSamplingShowerTuning::ConstantResolution(), G4cout, G4endl, GFlashSamplingShowerTuning::NoiseResolution(), GVFlashHomoShowerTuning::ParAveA1(), GVFlashHomoShowerTuning::ParAveA2(), GVFlashHomoShowerTuning::ParAveA3(), GVFlashHomoShowerTuning::ParAveT1(), GVFlashHomoShowerTuning::ParRC1(), GVFlashHomoShowerTuning::ParRC2(), GVFlashHomoShowerTuning::ParRC3(), GVFlashHomoShowerTuning::ParRC4(), GVFlashHomoShowerTuning::ParRho1(), GVFlashHomoShowerTuning::ParRho2(), GVFlashHomoShowerTuning::ParRT1(), GVFlashHomoShowerTuning::ParRT2(), GVFlashHomoShowerTuning::ParRT3(), GVFlashHomoShowerTuning::ParRT4(), GVFlashHomoShowerTuning::ParRT5(), GVFlashHomoShowerTuning::ParRT6(), GFlashSamplingShowerTuning::ParsAveA1(), GFlashSamplingShowerTuning::ParsAveT1(), GFlashSamplingShowerTuning::ParsAveT2(), GVFlashHomoShowerTuning::ParSigLogA1(), GVFlashHomoShowerTuning::ParSigLogA2(), GVFlashHomoShowerTuning::ParSigLogT1(), GVFlashHomoShowerTuning::ParSigLogT2(), GVFlashHomoShowerTuning::ParSpotA1(), GVFlashHomoShowerTuning::ParSpotA2(), GVFlashHomoShowerTuning::ParSpotN1(), GVFlashHomoShowerTuning::ParSpotN2(), GVFlashHomoShowerTuning::ParSpotT1(), GVFlashHomoShowerTuning::ParSpotT2(), GFlashSamplingShowerTuning::ParsRC1(), GFlashSamplingShowerTuning::ParsRC2(), GFlashSamplingShowerTuning::ParsRT1(), GFlashSamplingShowerTuning::ParsRT2(), GFlashSamplingShowerTuning::ParsWC1(), GFlashSamplingShowerTuning::ParsWC2(), GVFlashHomoShowerTuning::ParWC1(), GVFlashHomoShowerTuning::ParWC2(), GVFlashHomoShowerTuning::ParWC3(), GVFlashHomoShowerTuning::ParWC4(), GVFlashHomoShowerTuning::ParWC5(), GVFlashHomoShowerTuning::ParWC6(), GFlashSamplingShowerTuning::SamplingResolution(), and SetMaterial().

52  ParAveT2(0.), ParSigLogT1(0.), ParSigLogT2(0.),
53  ParSigLogA1(0.), ParSigLogA2(0.), ParRho1(0.), ParRho2(0.), ParsAveA2(0.),
54  AveLogAlphah(0.), AveLogTmaxh(0.), SigmaLogAlphah(0.), SigmaLogTmaxh(0.),
55  Rhoh(0.), Alphah(0.), Tmaxh(0.), Betah(0.), AveLogAlpha(0.), AveLogTmax(0.),
56  SigmaLogAlpha(0.), SigmaLogTmax(0.), Rho(0.), Alpha(0.), Tmax(0.), Beta(0.)
57 {
58  if(!aPar) { thePar = new GFlashSamplingShowerTuning; owning = true; }
59  else { thePar = aPar; owning = false; }
60 
61  SetMaterial(aMat1,aMat2 );
62  d1=dd1;
63  d2=dd2;
64 
65  // Longitudinal Coefficients for a homogenious calo
66 
67  // shower max
68  ParAveT1 = thePar->ParAveT1(); // ln (ln y -0.812)
69  ParAveA1 = thePar->ParAveA1(); // ln a (0.81 + (0.458 + 2.26/Z)ln y)
70  ParAveA2 = thePar->ParAveA2();
71  ParAveA3 = thePar->ParAveA3();
72  // Sampling
73  ParsAveT1 = thePar->ParsAveT1(); // T_sam = log(exp( log T_hom) + t1*Fs-1 + t2*(1-ehat));
74  ParsAveT2 = thePar->ParsAveT2();
75  ParsAveA1 = thePar->ParsAveA1();
76  // Variance of shower max sampling
77  ParsSigLogT1 = thePar->ParSigLogT1(); // Sigma T1 (-2.5 + 1.25 ln y)**-1
78  ParsSigLogT2 = thePar->ParSigLogT2();
79  // variance of 'alpha'
80  ParsSigLogA1 = thePar->ParSigLogA1(); // Sigma a (-0.82 + 0.79 ln y)**-1
81  ParsSigLogA2 = thePar->ParSigLogA2();
82  // correlation alpha%T
83  ParsRho1 = thePar->ParRho1(); // Rho = 0.784 -0.023 ln y
84  ParsRho2 = thePar->ParRho2();
85 
86  // Radial Coefficients
87  // r_C (tau)= z_1 +z_2 tau
88  // r_t (tau)= k1 (std::exp (k3(tau -k2 ))+std::exp (k_4 (tau- k_2))))
89  ParRC1 = thePar->ParRC1(); // z_1 = 0.0251 + 0.00319 ln E
90  ParRC2 = thePar->ParRC2();
91  ParRC3 = thePar->ParRC3(); // z_2 = 0.1162 + - 0.000381 Z
92  ParRC4 = thePar->ParRC4();
93 
94  ParWC1 = thePar->ParWC1();
95  ParWC2 = thePar->ParWC2();
96  ParWC3 = thePar->ParWC3();
97  ParWC4 = thePar->ParWC4();
98  ParWC5 = thePar->ParWC5();
99  ParWC6 = thePar->ParWC6();
100  ParRT1 = thePar->ParRT1();
101  ParRT2 = thePar->ParRT2();
102  ParRT3 = thePar->ParRT3();
103  ParRT4 = thePar->ParRT4();
104  ParRT5 = thePar->ParRT5();
105  ParRT6 = thePar->ParRT6();
106 
107  //additional sampling parameter
108  ParsRC1= thePar->ParsRC1();
109  ParsRC2= thePar->ParsRC2();
110  ParsWC1= thePar->ParsWC1();
111  ParsWC2= thePar->ParsWC2();
112  ParsRT1= thePar->ParsRT1();
113  ParsRT2= thePar->ParsRT2();
114 
115  // Coeff for fluctuedted radial profiles for a sampling media
116  ParsSpotT1 = thePar->ParSpotT1(); // T_spot = T_hom =(0.698 + 0.00212)
117  ParsSpotT2 = thePar->ParSpotT2();
118  ParsSpotA1 = thePar->ParSpotA1(); // a_spot= a_hom (0.639 + 0.00334)
119  ParsSpotA2 = thePar->ParSpotA2();
120  ParsSpotN1 = thePar->ParSpotN1(); // N_Spot 93 * ln(Z) E ** 0.876
121  ParsSpotN2 = thePar->ParSpotN2();
122  SamplingResolution = thePar->SamplingResolution();
123  ConstantResolution = thePar->ConstantResolution();
124  NoiseResolution = thePar->NoiseResolution();
125 
126  // Inits
127  NSpot = 0.00;
128  AlphaNSpot = 0.00;
129  TNSpot = 0.00;
130  BetaNSpot = 0.00;
131  RadiusCore = 0.00;
132  WeightCore = 0.00;
133  RadiusTail = 0.00;
135 
136  G4cout << "/********************************************/ " << G4endl;
137  G4cout << " - GFlashSamplingShowerParameterisation::Constructor - " << G4endl;
138  G4cout << "/********************************************/ " << G4endl;
139 }
G4GLOB_DLL std::ostream G4cout
void SetMaterial(G4Material *mat1, G4Material *mat2)
#define G4endl
Definition: G4ios.hh:61
GFlashSamplingShowerParameterisation::~GFlashSamplingShowerParameterisation ( )

Definition at line 143 of file GFlashSamplingShowerParameterisation.cc.

144 {
145  if(owning) { delete thePar; }
146 }

Member Function Documentation

G4double GFlashSamplingShowerParameterisation::ApplySampling ( const G4double  DEne,
const G4double  Energy 
)

Definition at line 294 of file GFlashSamplingShowerParameterisation.cc.

References G4INCL::DeJongSpin::shoot().

295 {
296  G4double DEneFluctuated = DEne;
297  G4double Resolution = std::pow(SamplingResolution,2);
298 
299  // +pow(NoiseResolution,2)/ //@@@@@@@@ FIXME
300  // Energy*(1.*MeV)+
301  // pow(ConstantResolution,2)*
302  // Energy/(1.*MeV);
303 
304  if(Resolution >0.0 && DEne > 0.00)
305  {
306  G4float x1=DEne/Resolution;
307  G4float x2 = G4RandGamma::shoot(x1, 1.0)*Resolution;
308  DEneFluctuated=x2;
309  }
310  return DEneFluctuated;
311 }
ThreeVector shoot(const G4int Ap, const G4int Af)
float G4float
Definition: G4Types.hh:77
double G4double
Definition: G4Types.hh:76
void GFlashSamplingShowerParameterisation::ComputeRadialParameters ( G4double  y,
G4double  Tau 
)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 383 of file GFlashSamplingShowerParameterisation.cc.

References python.hepunit::GeV.

Referenced by GenerateRadius().

384 {
385  G4double z1 = ParRC1 + ParRC2* std::log(Energy/GeV); //ok
386  G4double z2 = ParRC3+ParRC4*Zeff; //ok
387  RadiusCore = z1 + z2 * Tau; //ok
388  G4double p1 = ParWC1+ParWC2*Zeff; //ok
389  G4double p2 = ParWC3+ParWC4*Zeff; //ok
390  G4double p3 = ParWC5+ParWC6*std::log(Energy/GeV); //ok
391  WeightCore = p1 * std::exp( (p2-Tau)/p3- std::exp( (p2-Tau) /p3) ); //ok
392 
393  G4double k1 = ParRT1+ParRT2*Zeff; // ok
394  G4double k2 = ParRT3; // ok
395  G4double k3 = ParRT4; // ok
396  G4double k4 = ParRT5+ParRT6* std::log(Energy/GeV); // ok
397 
398  RadiusTail = k1*(std::exp(k3*(Tau-k2))
399  + std::exp(k4*(Tau-k2)) ); //ok
400 
401  // sampling calorimeter
402 
403  RadiusCore = RadiusCore + ParsRC1*(1-ehat) + ParsRC2/Fs*std::exp(-Tau); //ok
404  WeightCore = WeightCore + (1-ehat)
405  * (ParsWC1+ParsWC2/Fs * std::exp(-std::pow((Tau-1.),2))); //ok
406  RadiusTail = RadiusTail + (1-ehat)* ParsRT1+ ParsRT2/Fs *std::exp(-Tau); //ok
407 }
double G4double
Definition: G4Types.hh:76
G4double GFlashSamplingShowerParameterisation::ComputeTau ( G4double  LongitudinalPosition)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 372 of file GFlashSamplingShowerParameterisation.cc.

Referenced by GenerateRadius().

373 {
374  G4double tau = LongitudinalPosition / Tmax/ X0eff //<t> = T* a /(a - 1)
375  * (Alpha-1.00) /Alpha
376  * std::exp(AveLogAlpha)/(std::exp(AveLogAlpha)-1.); //ok
377  return tau;
378 }
double G4double
Definition: G4Types.hh:76
void GFlashSamplingShowerParameterisation::ComputeZAX0EFFetc ( )

Definition at line 174 of file GFlashSamplingShowerParameterisation.cc.

References python.hepunit::cm, python.hepunit::cm2, python.hepunit::cm3, g(), G4cout, G4endl, python.hepunit::MeV, and python.hepunit::mm.

Referenced by GFlashSamplingShowerParameterisation().

175 {
176  G4cout << "/************ ComputeZAX0EFFetc ************/" << G4endl;
177  G4cout << " - GFlashSamplingShowerParameterisation::Material - " << G4endl;
178 
179  G4double Es = 21*MeV; //constant
180 
181  // material and geometry parameters for a sampling calorimeter
182  G4double denominator = (d1*density1 + d2*density2);
183  G4double W1 = (d1*density1) / denominator;
184  G4double W2 = (d2*density2)/denominator;
185  Zeff = ( W1*Z2 ) + (W2*Z1); //X0*Es/Ec;
186  Aeff = ( W1*A1 ) + (W2*A2);
187  X0eff =(1/ (( W1 / X01) +( W2 / X02)));
188  Rhoeff = ( (d1 *density1 ) + (d2 * density2 ))/G4double (d2 + d1 );
189  Rmeff = 1/ ((((W1*Ec1)/ X01) + ((W2* Ec2)/ X02) ) / Es ) ;
190  Eceff = X0eff *((W1*Ec1)/ X01 + (W2* Ec2)/ X02 );
191  Fs = X0eff/G4double ((d1/mm )+(d2/mm) );
192  ehat = (1. / (1+ 0.007*(Z1- Z2)));
193 
194  G4cout << "W1= " << W1 << G4endl;
195  G4cout << "W2= " << W2 << G4endl;
196  G4cout << "effective quantities Zeff = "<<Zeff<< G4endl;
197  G4cout << "effective quantities Aeff = "<<Aeff<< G4endl;
198  G4cout << "effective quantities Rhoeff = "<<Rhoeff/g *cm3<<" g/cm3" << G4endl;
199  G4cout << "effective quantities X0eff = "<<X0eff/cm <<" cm" << G4endl;
200 
201  X0eff = X0eff * Rhoeff;
202 
203  G4cout << "effective quantities X0eff = "<<X0eff/g*cm2 <<" g/cm2" << G4endl;
204  X0eff = X0eff /Rhoeff;
205  G4cout << "effective quantities RMeff = "<<Rmeff/cm<<" cm" << G4endl;
206  Rmeff = Rmeff* Rhoeff;
207  G4cout << "effective quantities RMeff = "<<Rmeff/g *cm2<<" g/cm2" << G4endl;
208  Rmeff = Rmeff/ Rhoeff;
209  G4cout << "effective quantities Eceff = "<<Eceff/MeV<< " MeV"<< G4endl;
210  G4cout << "effective quantities Fs = "<<Fs<<G4endl;
211  G4cout << "effective quantities ehat = "<<ehat<<G4endl;
212  G4cout << "/********************************************/ " <<G4endl;
213 }
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GFlashSamplingShowerParameterisation::GenerateExponential ( G4double  Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 412 of file GFlashSamplingShowerParameterisation.cc.

References G4INCL::DeJongSpin::shoot().

413 {
414  G4double ParExp1 = 9./7.*X0eff;
415  G4double random = -ParExp1*G4RandExponential::shoot() ;
416  return random;
417 }
ThreeVector shoot(const G4int Ap, const G4int Af)
double G4double
Definition: G4Types.hh:76
void GFlashSamplingShowerParameterisation::GenerateLongitudinalProfile ( G4double  Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 218 of file GFlashSamplingShowerParameterisation.cc.

References FatalException, and G4Exception().

219 {
220  if ((material1==0) || (material2 ==0))
221  {
222  G4Exception("GFlashSamplingShowerParameterisation::GenerateLongitudinalProfile()",
223  "InvalidSetup", FatalException, "No material initialized!");
224  }
225  G4double y = Energy/Eceff;
226  ComputeLongitudinalParameters(y);
227  GenerateEnergyProfile(y);
228  GenerateNSpotProfile(y);
229 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76
G4double GFlashSamplingShowerParameterisation::GeneratePhi ( )
G4double GFlashSamplingShowerParameterisation::GenerateRadius ( G4int  ispot,
G4double  Energy,
G4double  LongitudinalPosition 
)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 342 of file GFlashSamplingShowerParameterisation.cc.

References ComputeRadialParameters(), ComputeTau(), DBL_MAX, G4UniformRand, and G4INCL::Math::min().

343 {
344  if(ispot < 1)
345  {
346  // Determine lateral parameters in the middle of the step.
347  // They depend on energy & position along step
348  //
349  G4double Tau = ComputeTau(LongitudinalPosition);
350  ComputeRadialParameters(Energy,Tau);
351  }
352 
353  G4double Radius;
354  G4double Random1 = G4UniformRand();
355  G4double Random2 = G4UniformRand();
356  if(Random1 <WeightCore) //WeightCore = p < w_i
357  {
358  Radius = Rmeff * RadiusCore * std::sqrt( Random2/(1. - Random2) );
359  }
360  else
361  {
362  Radius = Rmeff * RadiusTail * std::sqrt( Random2/(1. - Random2) );
363  }
364  Radius = std::min(Radius,DBL_MAX);
365  return Radius;
366 }
G4double ComputeTau(G4double LongitudinalPosition)
#define G4UniformRand()
Definition: Randomize.hh:87
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double GFlashSamplingShowerParameterisation::GetAveR90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 74 of file GFlashSamplingShowerParameterisation.hh.

74 {return (1.5 * Rmeff);} //ok
G4double GFlashSamplingShowerParameterisation::GetAveR99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 73 of file GFlashSamplingShowerParameterisation.hh.

73 {return (3.5 * Rmeff);}
G4double GFlashSamplingShowerParameterisation::GetAveT90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 78 of file GFlashSamplingShowerParameterisation.hh.

78 {return (2.5* X0eff* std::exp( AveLogTmax));}
G4double GFlashSamplingShowerParameterisation::GetAveT99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 77 of file GFlashSamplingShowerParameterisation.hh.

77 {return (X0eff*AveLogTmax/(AveLogAlpha-1.00));}
G4double GFlashSamplingShowerParameterisation::GetAveTmx ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 76 of file GFlashSamplingShowerParameterisation.hh.

76 {return (X0eff*std::exp(AveLogTmax));}
G4double GFlashSamplingShowerParameterisation::GetEc ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 82 of file GFlashSamplingShowerParameterisation.hh.

82 {return Eceff;}
G4double GFlashSamplingShowerParameterisation::GetNspot ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 80 of file GFlashSamplingShowerParameterisation.hh.

80 {return NSpot;}
G4double GFlashSamplingShowerParameterisation::GetRm ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 83 of file GFlashSamplingShowerParameterisation.hh.

83 {return Rmeff;}
G4double GFlashSamplingShowerParameterisation::GetX0 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 81 of file GFlashSamplingShowerParameterisation.hh.

81 {return X0eff;}
G4double GFlashSamplingShowerParameterisation::IntegrateEneLongitudinal ( G4double  LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 316 of file GFlashSamplingShowerParameterisation.cc.

References GVFlashShowerParameterisation::gam().

317 {
318  G4double LongitudinalStepInX0 = LongitudinalStep / X0eff;
319  G4float x1= Betah*LongitudinalStepInX0;
320  G4float x2= Alphah;
321  float x3 = gam(x1,x2);
322  G4double DEne=x3;
323  return DEne;
324 }
float G4float
Definition: G4Types.hh:77
G4double gam(G4double x, G4double a) const
double G4double
Definition: G4Types.hh:76
G4double GFlashSamplingShowerParameterisation::IntegrateNspLongitudinal ( G4double  LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 329 of file GFlashSamplingShowerParameterisation.cc.

References GVFlashShowerParameterisation::gam().

330 {
331  G4double LongitudinalStepInX0 = LongitudinalStep / X0eff;
332  G4float x1 = BetaNSpot*LongitudinalStepInX0;
333  G4float x2 = AlphaNSpot;
334  G4float x3 = gam(x1,x2);
335  G4double DNsp = x3;
336  return DNsp;
337 }
float G4float
Definition: G4Types.hh:77
G4double gam(G4double x, G4double a) const
double G4double
Definition: G4Types.hh:76
void GFlashSamplingShowerParameterisation::SetMaterial ( G4Material mat1,
G4Material mat2 
)

Definition at line 151 of file GFlashSamplingShowerParameterisation.cc.

References G4Material::GetDensity(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4Material::GetRadlen(), and python.hepunit::MeV.

Referenced by GFlashSamplingShowerParameterisation().

152 {
153  G4double Es = 21*MeV;
154  material1= mat1;
155  Z1 = GetEffZ(material1);
156  A1 = GetEffA(material1);
157  density1 = material1->GetDensity();
158  X01 = material1->GetRadlen();
159  Ec1 = 2.66 * std::pow((X01 * Z1 / A1),1.1);
160  Rm1 = X01*Es/Ec1;
161 
162  material2= mat2;
163  Z2 = GetEffZ(material2);
164  A2 = GetEffA(material2);
165  density2 = material2->GetDensity();
166  X02 = material2->GetRadlen();
167  Ec2 = 2.66 * std::pow((X02 * Z2 / A2),1.1);
168  Rm2 = X02*Es/Ec2;
169  // PrintMaterial();
170 }
G4double GetEffZ(const G4Material *material)
G4double GetDensity() const
Definition: G4Material.hh:178
G4double GetRadlen() const
Definition: G4Material.hh:218
G4double GetEffA(const G4Material *material)
double G4double
Definition: G4Types.hh:76

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