Geant4 11.2.2
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.

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

◆ ~GFlashHomoShowerParameterisation()

GFlashHomoShowerParameterisation::~GFlashHomoShowerParameterisation ( )

Definition at line 154 of file GFlashHomoShowerParameterisation.cc.

155{
156 delete thePar;
157}

Member Function Documentation

◆ ComputeRadialParameters()

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

Implements GVFlashShowerParameterisation.

Definition at line 274 of file GFlashHomoShowerParameterisation.cc.

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

Referenced by GenerateRadius().

◆ ComputeTau()

G4double GFlashHomoShowerParameterisation::ComputeTau ( G4double LongitudinalPosition)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 265 of file GFlashHomoShowerParameterisation.cc.

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

Referenced by GenerateRadius().

◆ ComputeZAX0EFFetc()

void GFlashHomoShowerParameterisation::ComputeZAX0EFFetc ( )

◆ GenerateExponential()

G4double GFlashHomoShowerParameterisation::GenerateExponential ( G4double Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 296 of file GFlashHomoShowerParameterisation.cc.

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

◆ GenerateLongitudinalProfile()

void GFlashHomoShowerParameterisation::GenerateLongitudinalProfile ( G4double Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 159 of file GFlashHomoShowerParameterisation.cc.

161{
162 if (material==0)
163 {
164 G4Exception("GFlashHomoShowerParameterisation::GenerateLongitudinalProfile()",
165 "InvalidSetup", FatalException, "No material initialized!");
166 }
167
168 G4double y = Energy/Ec;
169 ComputeLongitudinalParameters(y);
170 GenerateEnergyProfile(y);
171 GenerateNSpotProfile(y);
172}
@ 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 237 of file GFlashHomoShowerParameterisation.cc.

239{
240 if(ispot < 1)
241 {
242 // Determine lateral parameters in the middle of the step.
243 // They depend on energy & position along step.
244 //
245 G4double Tau = ComputeTau(LongitudinalPosition);
246 ComputeRadialParameters(Energy,Tau);
247 }
248
249 G4double Radius;
250 G4double Random1 = G4UniformRand();
251 G4double Random2 = G4UniformRand();
252
253 if(Random1 <WeightCore) //WeightCore = p < w_i
254 {
255 Radius = Rm * RadiusCore * std::sqrt( Random2/(1. - Random2) );
256 }
257 else
258 {
259 Radius = Rm * RadiusTail * std::sqrt( Random2/(1. - Random2) );
260 }
261 Radius = std::min(Radius,DBL_MAX);
262 return Radius;
263}
#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 72 of file GFlashHomoShowerParameterisation.hh.

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

◆ GetAveR99()

G4double GFlashHomoShowerParameterisation::GetAveR99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 71 of file GFlashHomoShowerParameterisation.hh.

71{return (3.5 * Rm);}

◆ GetAveT90()

G4double GFlashHomoShowerParameterisation::GetAveT90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 76 of file GFlashHomoShowerParameterisation.hh.

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

◆ GetAveT99()

G4double GFlashHomoShowerParameterisation::GetAveT99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 75 of file GFlashHomoShowerParameterisation.hh.

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

◆ GetAveTmx()

G4double GFlashHomoShowerParameterisation::GetAveTmx ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 74 of file GFlashHomoShowerParameterisation.hh.

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

◆ GetEc()

G4double GFlashHomoShowerParameterisation::GetEc ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 80 of file GFlashHomoShowerParameterisation.hh.

80{return Ec;}

◆ GetNspot()

G4double GFlashHomoShowerParameterisation::GetNspot ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 78 of file GFlashHomoShowerParameterisation.hh.

78{ return NSpot;}

◆ GetRm()

G4double GFlashHomoShowerParameterisation::GetRm ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 81 of file GFlashHomoShowerParameterisation.hh.

81{return Rm;}

◆ GetX0()

G4double GFlashHomoShowerParameterisation::GetX0 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 79 of file GFlashHomoShowerParameterisation.hh.

79{return X0;}

◆ IntegrateEneLongitudinal()

G4double GFlashHomoShowerParameterisation::IntegrateEneLongitudinal ( G4double LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 214 of file GFlashHomoShowerParameterisation.cc.

216{
217 G4double LongitudinalStepInX0 = LongitudinalStep / X0;
218 G4float x1= Betah*LongitudinalStepInX0;
219 G4float x2= Alphah;
220 float x3 = gam(x1,x2);
221 G4double DEne=x3;
222 return DEne;
223}
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 225 of file GFlashHomoShowerParameterisation.cc.

227{
228 G4double LongitudinalStepInX0 = LongitudinalStep / X0;
229 G4float x1 = BetaNSpot*LongitudinalStepInX0;
230 G4float x2 = AlphaNSpot;
231 G4float x3 = gam(x1,x2);
232 G4double DNsp = x3;
233 return DNsp;
234}

◆ SetMaterial()

void GFlashHomoShowerParameterisation::SetMaterial ( G4Material * mat)

Definition at line 141 of file GFlashHomoShowerParameterisation.cc.

142{
143 material= mat;
144 Z = GetEffZ(material);
145 A = GetEffA(material);
146 density = material->GetDensity()/(g/cm3);
147 X0 = material->GetRadlen();
148 Ec = 2.66 * std::pow((X0 * Z / A),1.1);
149 G4double Es = 21*MeV;
150 Rm = X0*Es/Ec;
151 // PrintMaterial();
152}
G4double GetDensity() const
G4double GetRadlen() const
G4double GetEffZ(const G4Material *material)
G4double GetEffA(const G4Material *material)

Referenced by GFlashHomoShowerParameterisation().


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