Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GFlashHomoShowerParameterisation Class Reference

#include <GFlashHomoShowerParameterisation.hh>

Inheritance diagram for GFlashHomoShowerParameterisation:
Collaboration 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 50 of file GFlashHomoShowerParameterisation.hh.

Constructor & Destructor Documentation

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

Definition at line 49 of file GFlashHomoShowerParameterisation.cc.

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

Here is the call graph for this function:

GFlashHomoShowerParameterisation::~GFlashHomoShowerParameterisation ( )

Definition at line 152 of file GFlashHomoShowerParameterisation.cc.

153 {
154  if(owning) { delete thePar; }
155 }

Member Function Documentation

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

Implements GVFlashShowerParameterisation.

Definition at line 273 of file GFlashHomoShowerParameterisation.cc.

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

Here is the caller graph for this function:

G4double GFlashHomoShowerParameterisation::ComputeTau ( G4double  LongitudinalPosition)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 264 of file GFlashHomoShowerParameterisation.cc.

265 {
266  G4double tau = LongitudinalPosition / Tmaxh / X0 //<t> = T* a /(a - 1)
267  * (Alphah-1.00) /Alphah *
268  std::exp(AveLogAlphah)/(std::exp(AveLogAlphah)-1.); //ok
269  return tau;
270 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

void GFlashHomoShowerParameterisation::ComputeZAX0EFFetc ( )
G4double GFlashHomoShowerParameterisation::GenerateExponential ( G4double  Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 295 of file GFlashHomoShowerParameterisation.cc.

296 {
297  G4double ParExp1 = 9./7.*X0;
298  G4double random = -ParExp1*G4RandExponential::shoot() ;
299  return random;
300 }
ThreeVector shoot(const G4int Ap, const G4int Af)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void GFlashHomoShowerParameterisation::GenerateLongitudinalProfile ( G4double  Energy)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 158 of file GFlashHomoShowerParameterisation.cc.

159 {
160  if (material==0)
161  {
162  G4Exception("GFlashHomoShowerParameterisation::GenerateLongitudinalProfile()",
163  "InvalidSetup", FatalException, "No material initialized!");
164  }
165 
166  G4double y = Energy/Ec;
167  ComputeLongitudinalParameters(y);
168  GenerateEnergyProfile(y);
169  GenerateNSpotProfile(y);
170 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements GVFlashShowerParameterisation.

Definition at line 236 of file GFlashHomoShowerParameterisation.cc.

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

Here is the call graph for this function:

G4double GFlashHomoShowerParameterisation::GetAveR90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 73 of file GFlashHomoShowerParameterisation.hh.

73 {return (1.5 * Rm);} //ok
G4double GFlashHomoShowerParameterisation::GetAveR99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 72 of file GFlashHomoShowerParameterisation.hh.

72 {return (3.5 * Rm);}
G4double GFlashHomoShowerParameterisation::GetAveT90 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 77 of file GFlashHomoShowerParameterisation.hh.

77 {return (2.5* X0*std::exp( AveLogTmaxh) );}
G4double GFlashHomoShowerParameterisation::GetAveT99 ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 76 of file GFlashHomoShowerParameterisation.hh.

76 {return (X0 * AveLogTmaxh/(AveLogAlphah-1.00));}
G4double GFlashHomoShowerParameterisation::GetAveTmx ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 75 of file GFlashHomoShowerParameterisation.hh.

75 {return (X0 * std::exp(AveLogTmaxh));}
G4double GFlashHomoShowerParameterisation::GetEc ( )
inlinevirtual
G4double GFlashHomoShowerParameterisation::GetNspot ( )
inlinevirtual

Implements GVFlashShowerParameterisation.

Definition at line 79 of file GFlashHomoShowerParameterisation.hh.

79 { return NSpot;}
G4double GFlashHomoShowerParameterisation::GetRm ( )
inlinevirtual
G4double GFlashHomoShowerParameterisation::GetX0 ( )
inlinevirtual
G4double GFlashHomoShowerParameterisation::IntegrateEneLongitudinal ( G4double  LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 213 of file GFlashHomoShowerParameterisation.cc.

214 {
215  G4double LongitudinalStepInX0 = LongitudinalStep / X0;
216  G4float x1= Betah*LongitudinalStepInX0;
217  G4float x2= Alphah;
218  float x3 = gam(x1,x2);
219  G4double DEne=x3;
220  return DEne;
221 }
float G4float
Definition: G4Types.hh:77
G4double gam(G4double x, G4double a) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double GFlashHomoShowerParameterisation::IntegrateNspLongitudinal ( G4double  LongitudinalStep)
virtual

Implements GVFlashShowerParameterisation.

Definition at line 224 of file GFlashHomoShowerParameterisation.cc.

225 {
226  G4double LongitudinalStepInX0 = LongitudinalStep / X0;
227  G4float x1 = BetaNSpot*LongitudinalStepInX0;
228  G4float x2 = AlphaNSpot;
229  G4float x3 = gam(x1,x2);
230  G4double DNsp = x3;
231  return DNsp;
232 }
float G4float
Definition: G4Types.hh:77
G4double gam(G4double x, G4double a) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void GFlashHomoShowerParameterisation::SetMaterial ( G4Material mat)

Definition at line 139 of file GFlashHomoShowerParameterisation.cc.

140 {
141  material= mat;
142  Z = GetEffZ(material);
143  A = GetEffA(material);
144  density = material->GetDensity()/(g/cm3);
145  X0 = material->GetRadlen();
146  Ec = 2.66 * std::pow((X0 * Z / A),1.1);
147  G4double Es = 21*MeV;
148  Rm = X0*Es/Ec;
149  // PrintMaterial();
150 }
G4double GetEffZ(const G4Material *material)
G4double GetDensity() const
Definition: G4Material.hh:180
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4double GetRadlen() const
Definition: G4Material.hh:220
static constexpr double cm3
Definition: G4SIunits.hh:121
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetEffA(const G4Material *material)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:


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