Geant4-11
GFlashHomoShowerParameterisation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// ------------------------------------------------------------
29// GEANT 4 class implementation
30//
31// ------- GFlashHomoShowerParameterisation -------
32//
33// Authors: E.Barberio & Joanna Weng - 9.11.2004
34// ------------------------------------------------------------
35
36#include <cmath>
37
41#include "G4SystemOfUnits.hh"
42#include "Randomize.hh"
43#include "G4ios.hh"
44#include "G4Material.hh"
45#include "G4MaterialTable.hh"
46
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) {
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)
75
76 // Variance of shower max
77 ParSigLogT1 = thePar->ParSigLogT1(); // Sigma T1 (-1.4 + 1.26 ln y)**-1
79
80 // variance of 'alpha'
81 //
82 ParSigLogA1 = thePar->ParSigLogA1(); // Sigma a (-0.58 + 0.86 ln y)**-1
84
85 // correlation alpha%T
86 //
87 ParRho1 = thePar->ParRho1(); // Rho = 0.705 -0.023 ln y
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)
118
119 ParSpotA1 = thePar->ParSpotA1(); // a_spot= a_hom (0.639 + 0.00334)
121
122 ParSpotN1 = thePar->ParSpotN1(); // N_Spot 93 * ln(Z) E ** 0.876
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}
140
142{
143 material= mat;
144 Z = GetEffZ(material);
145 A = GetEffA(material);
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}
153
155{
156 delete thePar;
157}
158
161{
162 if (material==0)
163 {
164 G4Exception("GFlashHomoShowerParameterisation::GenerateLongitudinalProfile()",
165 "InvalidSetup", FatalException, "No material initialized!");
166 }
167
168 G4double y = Energy/Ec;
172}
173
174void
176{
177 AveLogTmaxh = std::log(ParAveT1 + std::log(y));
178 //ok <ln T hom>
179 AveLogAlphah = std::log(ParAveA1 + (ParAveA2+ParAveA3/Z)*std::log(y));
180 //ok <ln alpha hom>
181
182 SigmaLogTmaxh = 1.00/( ParSigLogT1 + ParSigLogT2*std::log(y)) ;
183 //ok sigma (ln T hom)
184 SigmaLogAlphah = 1.00/( ParSigLogA1 + ParSigLogA2*std::log(y));
185 //ok sigma (ln alpha hom)
186 Rhoh = ParRho1+ParRho2*std::log(y); //ok
187}
188
190{
191 G4double Correlation1h = std::sqrt((1+Rhoh)/2);
192 G4double Correlation2h = std::sqrt((1-Rhoh)/2);
193
194 G4double Random1 = G4RandGauss::shoot();
195 G4double Random2 = G4RandGauss::shoot();
196
197 // Parameters for Enenrgy Profile including correaltion and sigmas
198
199 Tmaxh = std::exp( AveLogTmaxh + SigmaLogTmaxh *
200 (Correlation1h*Random1 + Correlation2h*Random2) );
201 Alphah = std::exp( AveLogAlphah + SigmaLogAlphah *
202 (Correlation1h*Random1 - Correlation2h*Random2) );
203 Betah = (Alphah-1.00)/Tmaxh;
204}
205
207{
208 TNSpot = Tmaxh * (ParSpotT1+ParSpotT2*Z); // ok
210 BetaNSpot = (AlphaNSpot-1.00)/TNSpot; // ok
211 NSpot = ParSpotN1 * std::log(Z)*std::pow((y*Ec)/GeV,ParSpotN2 ); // ok
212}
213
215IntegrateEneLongitudinal(G4double LongitudinalStep)
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}
224
226IntegrateNspLongitudinal(G4double LongitudinalStep)
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}
235
236
238GenerateRadius(G4int ispot, G4double Energy, G4double LongitudinalPosition)
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}
264
266ComputeTau(G4double LongitudinalPosition)
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}
273
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}
295
297GenerateExponential(const G4double /* Energy */ )
298{
299 G4double ParExp1 = 9./7.*X0;
300 G4double random = -ParExp1*G4RandExponential::shoot() ;
301 return random;
302}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double cm3
Definition: G4SIunits.hh:101
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
float G4float
Definition: G4Types.hh:84
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetDensity() const
Definition: G4Material.hh:176
G4double GetRadlen() const
Definition: G4Material.hh:216
G4double ComputeTau(G4double LongitudinalPosition)
G4double IntegrateEneLongitudinal(G4double LongitudinalStep)
GFlashHomoShowerParameterisation(G4Material *aMat, GVFlashHomoShowerTuning *aPar=0)
G4double GenerateRadius(G4int ispot, G4double Energy, G4double LongitudinalPosition)
G4double IntegrateNspLongitudinal(G4double LongitudinalStep)
void ComputeRadialParameters(G4double y, G4double Tau)
G4double GetEffZ(const G4Material *material)
G4double GetEffA(const G4Material *material)
G4double gam(G4double x, G4double a) const
ThreeVector shoot(const G4int Ap, const G4int Af)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define DBL_MAX
Definition: templates.hh:62