Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
GFlashHomoShowerParameterisation Class Reference

#include <GFlashHomoShowerParameterisation.hh>

+ Inheritance diagram for GFlashHomoShowerParameterisation:

Public Member Functions

 GFlashHomoShowerParameterisation (G4Material *aMat, GVFlashHomoShowerTuning *aPar=0)
 
 ~GFlashHomoShowerParameterisation ()
 
void ComputeRadialParameters (G4double y, G4double Tau)
 
void GenerateLongitudinalProfile (G4double Energy)
 
void ComputeZAX0EFFetc ()
 
G4double IntegrateEneLongitudinal (G4double LongitudinalStep)
 
G4double IntegrateNspLongitudinal (G4double LongitudinalStep)
 
G4double ComputeTau (G4double LongitudinalPosition)
 
G4double GeneratePhi ()
 
G4double GenerateRadius (G4int ispot, G4double Energy, G4double LongitudinalPosition)
 
G4double GenerateExponential (G4double Energy)
 
void SetMaterial (G4Material *mat)
 
G4double GetAveR99 ()
 
G4double GetAveR90 ()
 
G4double GetAveTmx ()
 
G4double GetAveT99 ()
 
G4double GetAveT90 ()
 
G4double GetNspot ()
 
G4double GetX0 ()
 
G4double GetEc ()
 
G4double GetRm ()
 
- 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 49 of file GFlashHomoShowerParameterisation.hh.

Constructor & Destructor Documentation

◆ GFlashHomoShowerParameterisation()

GFlashHomoShowerParameterisation::GFlashHomoShowerParameterisation ( G4Material * aMat,
GVFlashHomoShowerTuning * aPar = 0 )

Definition at line 47 of file GFlashHomoShowerParameterisation.cc.

50 ConstantResolution(0.),
51 NoiseResolution(0.),
52 SamplingResolution(0.),
53 AveLogAlphah(0.),
54 AveLogTmaxh(0.),
55 SigmaLogAlphah(0.),
56 SigmaLogTmaxh(0.),
57 Rhoh(0.),
58 Alphah(0.),
59 Tmaxh(0.),
60 Betah(0.)
61
62{
63 if (!aPar) {
64 thePar = new GVFlashHomoShowerTuning;
65 }
66 else {
67 thePar = aPar;
68 }
69
70 SetMaterial(aMat);
71 PrintMaterial(aMat);
72
73 /********************************************/
74 /* Homo Calorimeter */
75 /********************************************/
76 // Longitudinal Coefficients for a homogenious calo
77 // shower max
78 //
79 ParAveT1 = thePar->ParAveT1(); // ln (ln y -0.812)
80 ParAveA1 = thePar->ParAveA1(); // ln a (0.81 + (0.458 + 2.26/Z)ln y)
81 ParAveA2 = thePar->ParAveA2();
82 ParAveA3 = thePar->ParAveA3();
83
84 // Variance of shower max
85 ParSigLogT1 = thePar->ParSigLogT1(); // Sigma T1 (-1.4 + 1.26 ln y)**-1
86 ParSigLogT2 = thePar->ParSigLogT2();
87
88 // variance of 'alpha'
89 //
90 ParSigLogA1 = thePar->ParSigLogA1(); // Sigma a (-0.58 + 0.86 ln y)**-1
91 ParSigLogA2 = thePar->ParSigLogA2();
92
93 // correlation alpha%T
94 //
95 ParRho1 = thePar->ParRho1(); // Rho = 0.705 -0.023 ln y
96 ParRho2 = thePar->ParRho2();
97
98 // Radial Coefficients
99 // r_C (tau)= z_1 +z_2 tau
100 // r_t (tau)= k1 (std::exp (k3(tau -k2 ))+std::exp (k_4 (tau- k_2))))
101 //
102 ParRC1 = thePar->ParRC1(); // z_1 = 0.0251 + 0.00319 ln E
103 ParRC2 = thePar->ParRC2();
104
105 ParRC3 = thePar->ParRC3(); // z_2 = 0.1162 + - 0.000381 Z
106 ParRC4 = thePar->ParRC4();
107
108 ParWC1 = thePar->ParWC1();
109 ParWC2 = thePar->ParWC2();
110 ParWC3 = thePar->ParWC3();
111 ParWC4 = thePar->ParWC4();
112 ParWC5 = thePar->ParWC5();
113 ParWC6 = thePar->ParWC6();
114
115 ParRT1 = thePar->ParRT1();
116 ParRT2 = thePar->ParRT2();
117 ParRT3 = thePar->ParRT3();
118 ParRT4 = thePar->ParRT4();
119 ParRT5 = thePar->ParRT5();
120 ParRT6 = thePar->ParRT6();
121
122 // Coeff for fluctueted radial profiles for a uniform media
123 //
124 ParSpotT1 = thePar->ParSpotT1(); // T_spot = T_hom =(0.698 + 0.00212)
125 ParSpotT2 = thePar->ParSpotT2();
126
127 ParSpotA1 = thePar->ParSpotA1(); // a_spot= a_hom (0.639 + 0.00334)
128 ParSpotA2 = thePar->ParSpotA2();
129
130 ParSpotN1 = thePar->ParSpotN1(); // N_Spot 93 * ln(Z) E ** 0.876
131 ParSpotN2 = thePar->ParSpotN2();
132
133 // Inits
134
135 NSpot = 0.00;
136 AlphaNSpot = 0.00;
137 TNSpot = 0.00;
138 BetaNSpot = 0.00;
139
140 RadiusCore = 0.00;
141 WeightCore = 0.00;
142 RadiusTail = 0.00;
143
144 G4cout << "/********************************************/ " << G4endl;
145 G4cout << " - GFlashHomoShowerParameterisation::Constructor - " << G4endl;
146 G4cout << "/********************************************/ " << G4endl;
147}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ ~GFlashHomoShowerParameterisation()

GFlashHomoShowerParameterisation::~GFlashHomoShowerParameterisation ( )

Definition at line 165 of file GFlashHomoShowerParameterisation.cc.

166{
167 delete thePar;
168}

Member Function Documentation

◆ ComputeRadialParameters()

void GFlashHomoShowerParameterisation::ComputeRadialParameters ( G4double y,
G4double Tau )
virtual

Implements GVFlashShowerParameterisation.

Definition at line 277 of file GFlashHomoShowerParameterisation.cc.

278{
279 G4double z1 = ParRC1 + ParRC2 * std::log(Energy / GeV); // ok
280 G4double z2 = ParRC3 + ParRC4 * Z; // ok
281 RadiusCore = z1 + z2 * Tau; // ok
282
283 G4double p1 = ParWC1 + ParWC2 * Z; // ok
284 G4double p2 = ParWC3 + ParWC4 * Z; // ok
285 G4double p3 = ParWC5 + ParWC6 * std::log(Energy / GeV); // ok
286
287 WeightCore = p1 * std::exp((p2 - Tau) / p3 - std::exp((p2 - Tau) / p3)); // ok
288
289 G4double k1 = ParRT1 + ParRT2 * Z; // ok
290 G4double k2 = ParRT3; // ok
291 G4double k3 = ParRT4; // ok
292 G4double k4 = ParRT5 + ParRT6 * std::log(Energy / GeV); // ok
293
294 RadiusTail = k1 * (std::exp(k3 * (Tau - k2)) + std::exp(k4 * (Tau - k2))); // ok
295}
double G4double
Definition G4Types.hh:83

Referenced by GenerateRadius().

◆ ComputeTau()

G4double GFlashHomoShowerParameterisation::ComputeTau ( G4double LongitudinalPosition)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 269 of file GFlashHomoShowerParameterisation.cc.

270{
271 G4double tau = LongitudinalPosition / Tmaxh / X0 //<t> = T* a /(a - 1)
272 * (Alphah - 1.00) / Alphah * std::exp(AveLogAlphah)
273 / (std::exp(AveLogAlphah) - 1.); // ok
274 return tau;
275}

Referenced by GenerateRadius().

◆ ComputeZAX0EFFetc()

void GFlashHomoShowerParameterisation::ComputeZAX0EFFetc ( )

◆ GenerateExponential()

G4double GFlashHomoShowerParameterisation::GenerateExponential ( G4double Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 297 of file GFlashHomoShowerParameterisation.cc.

299{
300 G4double ParExp1 = 9./7.*X0;
301 G4double random = -ParExp1*G4RandExponential::shoot() ;
302 return random;
303}

◆ GenerateLongitudinalProfile()

void GFlashHomoShowerParameterisation::GenerateLongitudinalProfile ( G4double Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 170 of file GFlashHomoShowerParameterisation.cc.

171{
172 if (material == 0) {
173 G4Exception("GFlashHomoShowerParameterisation::GenerateLongitudinalProfile()", "InvalidSetup",
174 FatalException, "No material initialized!");
175 }
176
177 G4double y = Energy / Ec;
178 ComputeLongitudinalParameters(y);
179 GenerateEnergyProfile(y);
180 GenerateNSpotProfile(y);
181}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

◆ GeneratePhi()

G4double GFlashHomoShowerParameterisation::GeneratePhi ( )

◆ GenerateRadius()

G4double GFlashHomoShowerParameterisation::GenerateRadius ( G4int ispot,
G4double Energy,
G4double LongitudinalPosition )
virtual

Implements GVFlashShowerParameterisation.

Definition at line 243 of file GFlashHomoShowerParameterisation.cc.

245{
246 if (ispot < 1) {
247 // Determine lateral parameters in the middle of the step.
248 // They depend on energy & position along step.
249 //
250 G4double Tau = ComputeTau(LongitudinalPosition);
251 ComputeRadialParameters(Energy, Tau);
252 }
253
254 G4double Radius;
255 G4double Random1 = G4UniformRand();
256 G4double Random2 = G4UniformRand();
257
258 if (Random1 < WeightCore) // WeightCore = p < w_i
259 {
260 Radius = Rm * RadiusCore * std::sqrt(Random2 / (1. - Random2));
261 }
262 else {
263 Radius = Rm * RadiusTail * std::sqrt(Random2 / (1. - Random2));
264 }
265 Radius = std::min(Radius, DBL_MAX);
266 return Radius;
267}
#define G4UniformRand()
Definition Randomize.hh:52
G4double ComputeTau(G4double LongitudinalPosition)
void ComputeRadialParameters(G4double y, G4double Tau)
#define DBL_MAX
Definition templates.hh:62

◆ GetAveR90()

G4double GFlashHomoShowerParameterisation::GetAveR90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 69 of file GFlashHomoShowerParameterisation.hh.

69{ return (1.5 * Rm); } // ok

◆ GetAveR99()

G4double GFlashHomoShowerParameterisation::GetAveR99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 68 of file GFlashHomoShowerParameterisation.hh.

68{ return (3.5 * Rm); }

◆ GetAveT90()

G4double GFlashHomoShowerParameterisation::GetAveT90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 73 of file GFlashHomoShowerParameterisation.hh.

73{return (2.5* X0*std::exp( AveLogTmaxh) );}

◆ GetAveT99()

G4double GFlashHomoShowerParameterisation::GetAveT99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 72 of file GFlashHomoShowerParameterisation.hh.

72{return (X0 * AveLogTmaxh/(AveLogAlphah-1.00));}

◆ GetAveTmx()

G4double GFlashHomoShowerParameterisation::GetAveTmx ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 71 of file GFlashHomoShowerParameterisation.hh.

71{return (X0 * std::exp(AveLogTmaxh));}

◆ GetEc()

G4double GFlashHomoShowerParameterisation::GetEc ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 77 of file GFlashHomoShowerParameterisation.hh.

77{ return Ec; }

◆ GetNspot()

G4double GFlashHomoShowerParameterisation::GetNspot ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 75 of file GFlashHomoShowerParameterisation.hh.

75{ return NSpot;}

◆ GetRm()

G4double GFlashHomoShowerParameterisation::GetRm ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 78 of file GFlashHomoShowerParameterisation.hh.

78{ return Rm; }

◆ GetX0()

G4double GFlashHomoShowerParameterisation::GetX0 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 76 of file GFlashHomoShowerParameterisation.hh.

76{ return X0; }

◆ IntegrateEneLongitudinal()

G4double GFlashHomoShowerParameterisation::IntegrateEneLongitudinal ( G4double LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 222 of file GFlashHomoShowerParameterisation.cc.

224{
225 G4double LongitudinalStepInX0 = LongitudinalStep / X0;
226 G4float x1= Betah*LongitudinalStepInX0;
227 G4float x2= Alphah;
228 float x3 = gam(x1,x2);
229 G4double DEne=x3;
230 return DEne;
231}
float G4float
Definition G4Types.hh:84
G4double gam(G4double x, G4double a) const

◆ IntegrateNspLongitudinal()

G4double GFlashHomoShowerParameterisation::IntegrateNspLongitudinal ( G4double LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 233 of file GFlashHomoShowerParameterisation.cc.

234{
235 G4double LongitudinalStepInX0 = LongitudinalStep / X0;
236 G4float x1 = BetaNSpot * LongitudinalStepInX0;
237 G4float x2 = AlphaNSpot;
238 G4float x3 = gam(x1, x2);
239 G4double DNsp = x3;
240 return DNsp;
241}

◆ SetMaterial()

void GFlashHomoShowerParameterisation::SetMaterial ( G4Material * mat)

Definition at line 149 of file GFlashHomoShowerParameterisation.cc.

150{
151 material = mat;
152 Z = GetEffZ(material);
153 A = GetEffA(material);
154 density = material->GetDensity() / (g / cm3);
155 X0 = material->GetRadlen();
156 // O. I. Dovzhenkko and A. A. Pommanskii
157 Ec = 2.66 * std::pow((X0 * Z / A), 1.1);
158 // // Rossi appriximation
159 // Ec = 610.0 * MeV / (Z + 1.24);
160 const G4double Es = 21.2 * MeV;
161 Rm = X0 * Es / Ec;
162 // PrintMaterial();
163}
G4double GetEffZ(const G4Material *material)
G4double GetEffA(const G4Material *material)

Referenced by GFlashHomoShowerParameterisation().


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