Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Par03EMShowerModel Class Reference

Example fast simulation model for EM showers. More...

#include <Doxymodules_parameterisations.h>

Inheritance diagram for Par03EMShowerModel:
G4VFastSimulationModel

Public Member Functions

 Par03EMShowerModel (G4String, G4Region *)
 
 Par03EMShowerModel (G4String)
 
 ~Par03EMShowerModel ()
 
virtual G4bool ModelTrigger (const G4FastTrack &) final
 There are no kinematics constraints. True is returned.
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) final
 Model is applicable to electrons, positrons, and photons.
 
virtual void DoIt (const G4FastTrack &, G4FastStep &) final
 Take particle out of the full simulation (kill it at the entrance depositing all the energy).
 
void Print () const
 Print current settings.
 
void SetSigma (const G4double aSigma)
 Set standard deviation of a Gaussian distribution that describes the transverse shower profile.
 
G4double GetSigma () const
 Get standard deviation of a Gaussian distribution that describes the transverse shower profile.
 
void SetAlpha (const G4double aAlpha)
 Set alpha parameter of a Gamma distribution that describes the longitudinal shower profile.
 
G4double GetAlpha () const
 Get alpha parameter of a Gamma distribution that describes the longitudinal shower profile.
 
void SetBeta (const G4double aBeta)
 Set beta parameter of a Gamma distribution that describes the longitudinal shower profile.
 
G4double GetBeta () const
 Get beta parameter of a Gamma distribution that describes the longitudinal shower profile.
 
void SetNbOfHits (const G4int aNumber)
 Set number of (same energy) hits created in the parametrisation.
 
G4int GetNbOfHits () const
 Get number of (same energy) hits created in the parametrisation.s.
 
void SetLongMaxDepth (const G4double aDepth)
 Set maximum depth of shower created in fast simulation.
 
G4double GetLongMaxDepth () const
 Get maximum depth of shower created in fast simulation.
 

Private Member Functions

G4double Gamma (G4double x, G4double alpha, G4double beta)
 Gamma distribution.
 
G4double Gaussian (G4double x, G4double sigma=1, G4double x0=0)
 Gaussian distribution.
 

Private Attributes

Par03EMShowerMessengerfMessenger
 Messenger for configuration.
 
std::unique_ptr< G4FastSimHitMakerfHitMaker
 Helper class for creation of hits within the sensitive detector.
 
G4double fSigma = -1
 Standard deviation of the Gaussian distribution Can be changed with UI command /Par03/fastSim/transverseProfile/sigma <sigma> If sigma is smaller than 0, it will be estimated from the detector material (Moliere radius).
 
G4double fAlpha = -1
 Alpha parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudunalProfile/alpha <alpha> If alpha is smaller than 0, it will be estimated from particle energy and the detector material.
 
G4double fBeta = 0.5
 Beta parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudinalProfile/beta <beta>
 
G4int fNbOfHits = 100
 Number of (same energy) hits created by the parametrisation.
 
G4double fLongMaxDepth = 30
 Maximum depth of a shower created in fast simulation.
 

Detailed Description

Example fast simulation model for EM showers.

Parametrisation of electrons, positrons, and gammas. It is triggered if those particles enter the detector so that there is sufficient length for the shower development (max depth, controlled by the UI command).

Parametrisation is based on the PDG chapter on the electromagnetic cascades (chapter 33.5). Longitudinal profile of the shower is described with Gamma distribution, with beta parameter on average equal to 0.5 (default value, Fig. 33.21), and alpha parameter calcluated from the incident particle energy and material of the detector (critical energy) following Eq.(33.36).

Transverse profile is in this model approximated by the Gaussian distribution, with the mean along the shower axis (incident particle momentum direction) and the standard deviation calculated from the detector material (Moliere radius). This assumes that EM shower is in 90% contained within a cylinder of radius equal to Moliere radius, and that area below Gaussian distribution from mean-1.645 sigma to mean+1.645 sigma is also equal to 90% of total distribution.

Parameters of both distributions (alpha, beta for Gamma, sigma for Gaussian) can be overwritten by UI commands.

Parametrisation creates N hits of same energy (N can be set by UI command), using rejection sampling to generate position along shower axis from Gamma distribution, and then sampling from uniform and Gaussian distributions to sample phi and radius, respectively. Created hits are deposited in the detector using its readout geometry, using the helper class G4FastSimHitMaker that locates the volume, and calls appropriate sensitive detector class.

PDG Chapter 33: https://pdg.lbl.gov/2019/reviews/rpp2018-rev-passage-particles-matter.pdf

Definition at line 63 of file Doxymodules_parameterisations.h.

Constructor & Destructor Documentation

◆ Par03EMShowerModel() [1/2]

Par03EMShowerModel::Par03EMShowerModel ( G4String  aModelName,
G4Region aEnvelope 
)

Definition at line 38 of file Par03EMShowerModel.cc.

39 : G4VFastSimulationModel(aModelName, aEnvelope)
42{}
Messenger for the example fast simulation model.
Par03EMShowerMessenger * fMessenger
Messenger for configuration.
std::unique_ptr< G4FastSimHitMaker > fHitMaker
Helper class for creation of hits within the sensitive detector.

◆ Par03EMShowerModel() [2/2]

Par03EMShowerModel::Par03EMShowerModel ( G4String  aModelName)

Definition at line 46 of file Par03EMShowerModel.cc.

47 : G4VFastSimulationModel(aModelName)
50{}

◆ ~Par03EMShowerModel()

Par03EMShowerModel::~Par03EMShowerModel ( )
default

Member Function Documentation

◆ ModelTrigger()

G4bool Par03EMShowerModel::ModelTrigger ( const G4FastTrack aFastTrack)
finalvirtual

There are no kinematics constraints. True is returned.

Definition at line 68 of file Par03EMShowerModel.cc.

69{
70 // Check energy
71 if(aFastTrack.GetPrimaryTrack()->GetKineticEnergy() < 1 * GeV)
72 {
73 return false;
74 }
75 // Check length of detector
76 // Calculate depth of the detector along shower axis to verify if shower
77 // will fit inside. Required max shower depth is defined by fLongMaxDepth, and
78 // can be changed with UI command `/Par03/fastSim/longitudinalProfile/maxDepth
79 G4double X0 = aFastTrack.GetPrimaryTrack()->GetMaterial()->GetRadlen();
80 auto particleDirection = aFastTrack.GetPrimaryTrackLocalDirection();
81 auto particlePosition = aFastTrack.GetPrimaryTrackLocalPosition();
82 G4double detectorDepthInMM = aFastTrack.GetEnvelopeSolid()->DistanceToOut(
83 particlePosition, particleDirection);
84 G4double detectorDepthInX0 = detectorDepthInMM / X0;
85 // check if detector depth is sufficient to create showers
86 if(detectorDepthInX0 < fLongMaxDepth)
87 {
88 return false;
89 }
90 return true;
91}
G4double fLongMaxDepth
Maximum depth of a shower created in fast simulation.

◆ IsApplicable()

G4bool Par03EMShowerModel::IsApplicable ( const G4ParticleDefinition aParticleType)
finalvirtual

Model is applicable to electrons, positrons, and photons.

Definition at line 58 of file Par03EMShowerModel.cc.

60{
61 return &aParticleType == G4Electron::ElectronDefinition() ||
62 &aParticleType == G4Positron::PositronDefinition() ||
63 &aParticleType == G4Gamma::GammaDefinition();
64}

◆ DoIt()

void Par03EMShowerModel::DoIt ( const G4FastTrack aFastTrack,
G4FastStep aFastStep 
)
finalvirtual

Take particle out of the full simulation (kill it at the entrance depositing all the energy).

Calculate energy deposited in the detector according to Gamma distribution (along the particle direction) and Gaussian distribution in the transverse direction. Mean of the Gaussian is centred on the shower axis. Create energy deposits on a cylindrical mesh. Parameters of the mesh (size, number of cells) and of the distributions (alpha, beta for Gamma, sigma for Gaussian) can be set with UI commands.

Definition at line 95 of file Par03EMShowerModel.cc.

97{
98 // Remove particle from further processing by G4
99 aFastStep.KillPrimaryTrack();
100 aFastStep.ProposePrimaryTrackPathLength(0.0);
101 G4double energy = aFastTrack.GetPrimaryTrack()->GetKineticEnergy();
102 // No need to create any deposit, it will be handled by this model (and
103 // G4FastSimHitMaker that will call the sensitive detector)
104 aFastStep.ProposeTotalEnergyDeposited(0);
105 auto particlePosition = aFastTrack.GetPrimaryTrackLocalPosition();
106 auto particleDirection = aFastTrack.GetPrimaryTrackLocalDirection();
107
108 // Calculate how to create energy deposits
109 // Following PDG 33.5 chapter
110 // material calculation assumes homogeneous detector (true for Par03 example)
111 auto material = aFastTrack.GetPrimaryTrack()->GetMaterial();
112 G4double materialX0 = material->GetRadlen();
113 G4double materialZ = material->GetZ();
114 // EC estimation follows PDG fit to solids in Fig. 33.14 (rms 2.2%)
115 G4double materialEc = 610 * MeV / (materialZ + 1.24);
116 // RM estimation follows PDG Eq. (33.37) (rms 2.2%)
117 G4double materialRM = 21.2052 * MeV * materialX0 / materialEc;
118 G4double particleY = energy / materialEc;
119 // Estimate shower maximum and alpha parameter of Gamma distribution
120 // that describes the longitudinal profile (PDG Eq. (33.35))
121 // unless alpha is specified by UI command
122 if(fAlpha < 0)
123 {
124 // from PDG Eq. (33.36)
125 G4double particleTmax = std::log(particleY);
126 if(aFastTrack.GetPrimaryTrack()->GetParticleDefinition() ==
127 G4Gamma::GammaDefinition())
128 {
129 particleTmax += 0.5;
130 }
131 else
132 {
133 particleTmax -= 0.5;
134 }
135 fAlpha = particleTmax * fBeta + 1;
136 }
137 // Unless sigma of Gaussian distribution describing the transverse profile
138 // is specified by UI command, use value calculated from Moliere Radius
139 if(fSigma < 0)
140 {
141 // 90% of shower is contained within 1 * R_M
142 // 1.645 * std dev of Gaussian contains 90%
143 fSigma = materialRM / 1.645;
144 }
145
146 // Calculate rotation matrix along the particle momentum direction
147 // It will rotate the shower axes to match the incoming particle direction
148 G4RotationMatrix rotMatrix = G4RotationMatrix();
149 double particleTheta = particleDirection.theta();
150 double particlePhi = particleDirection.phi();
151 double epsilon = 1e-3;
152 rotMatrix.rotateY(particleTheta);
153 // do not use (random) phi if x==y==0
154 if(!(std::fabs(particleDirection.x()) < epsilon &&
155 std::fabs(particleDirection.y()) < epsilon))
156 rotMatrix.rotateZ(particlePhi);
157
158 // Create hits
159 // First use rejecton sampling to sample from Gamma distribution
160 // then get random numbers from uniform distribution for azimuthal angle, and
161 // from Gaussian for radius
162 G4ThreeVector position;
163 G4double gammaMax = Gamma((fAlpha - 1) / fBeta, fAlpha, fBeta);
164 G4int generatedHits = 0;
165 while(generatedHits < fNbOfHits)
166 {
167 G4double random1 = G4UniformRand() * fLongMaxDepth;
168 G4double random2 = G4UniformRand() * gammaMax;
169 if(Gamma(random1, fAlpha, fBeta) >= random2)
170 {
171 // Generate corresponding rho (phi) from Gaussian (flat) distribution
172 G4double phiPosition = G4UniformRand() * 2 * CLHEP::pi;
173 G4double rhoPosition = G4RandGauss::shoot(0, fSigma);
174 position = particlePosition +
175 rotMatrix * G4ThreeVector(rhoPosition * std::sin(phiPosition),
176 rhoPosition * std::cos(phiPosition),
177 random1 * materialX0);
178 // Create energy deposit in the detector
179 // This will call appropriate sensitive detector class
180 fHitMaker->make(G4FastHit(position, energy / fNbOfHits), aFastTrack);
181 generatedHits++;
182 }
183 }
184}
G4int fNbOfHits
Number of (same energy) hits created by the parametrisation.
G4double fBeta
Beta parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudinalPr...
G4double fAlpha
Alpha parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudunalP...
G4double fSigma
Standard deviation of the Gaussian distribution Can be changed with UI command /Par03/fastSim/transve...
G4double Gamma(G4double x, G4double alpha, G4double beta)
Gamma distribution.

◆ Print()

void Par03EMShowerModel::Print ( ) const

Print current settings.

Definition at line 188 of file Par03EMShowerModel.cc.

189{
190 G4cout << "Par03EMShowerModel: " << G4endl;
191 G4cout << "Gaussian distribution (transverse plane): \tmu = 0, sigma = "
192 << G4BestUnit(fSigma, "Length") << G4endl;
193 if(fSigma < 0)
194 G4cout << "Negative sigma value means that it will be recalculated "
195 "from the value of the Moliere radius of the detector material, "
196 "taking into account that 90% of the area below the Gaussian "
197 "distribution (from mu - 1.645 sigma to mu + 1.645 sigma) "
198 "corresponds to area within 1 Moliere radius."
199 << G4endl;
200 G4cout << "Gamma distribution (along shower axis): \talpha = " << fAlpha
201 << ", beta = " << fBeta << ", max depth = " << fLongMaxDepth << " X0"
202 << G4endl;
203 if(fAlpha < 0)
204 G4cout << "Negative alpha value means that it will be recalculated "
205 "from the critical energy of the detector material, particle "
206 "type, and beta parameter.\n alpha = beta * T_max, where T_max = "
207 "ln(E/E_C) + C\n where E is particle energy, E_C is critical "
208 "energy in the material, and constant C = -0.5 for electrons and "
209 "0.5 for photons (Eq. (33.36) from PDG)."
210 << G4endl;
211 G4cout << "Number of created energy deposits: " << fNbOfHits << G4endl;
212}

◆ SetSigma()

void Par03EMShowerModel::SetSigma ( const G4double  aSigma)
inline

Set standard deviation of a Gaussian distribution that describes the transverse shower profile.

Definition at line 94 of file Par03EMShowerModel.hh.

94{ fSigma = aSigma; };

◆ GetSigma()

G4double Par03EMShowerModel::GetSigma ( ) const
inline

Get standard deviation of a Gaussian distribution that describes the transverse shower profile.

Definition at line 97 of file Par03EMShowerModel.hh.

97{ return fSigma; };

◆ SetAlpha()

void Par03EMShowerModel::SetAlpha ( const G4double  aAlpha)
inline

Set alpha parameter of a Gamma distribution that describes the longitudinal shower profile.

Definition at line 100 of file Par03EMShowerModel.hh.

100{ fAlpha = aAlpha; };

◆ GetAlpha()

G4double Par03EMShowerModel::GetAlpha ( ) const
inline

Get alpha parameter of a Gamma distribution that describes the longitudinal shower profile.

Definition at line 103 of file Par03EMShowerModel.hh.

103{ return fAlpha; };

◆ SetBeta()

void Par03EMShowerModel::SetBeta ( const G4double  aBeta)
inline

Set beta parameter of a Gamma distribution that describes the longitudinal shower profile.

Definition at line 106 of file Par03EMShowerModel.hh.

106{ fBeta = aBeta; };

◆ GetBeta()

G4double Par03EMShowerModel::GetBeta ( ) const
inline

Get beta parameter of a Gamma distribution that describes the longitudinal shower profile.

Definition at line 109 of file Par03EMShowerModel.hh.

109{ return fBeta; };

◆ SetNbOfHits()

void Par03EMShowerModel::SetNbOfHits ( const G4int  aNumber)
inline

Set number of (same energy) hits created in the parametrisation.

Definition at line 111 of file Par03EMShowerModel.hh.

111{ fNbOfHits = aNumber; };

◆ GetNbOfHits()

G4int Par03EMShowerModel::GetNbOfHits ( ) const
inline

Get number of (same energy) hits created in the parametrisation.s.

Definition at line 113 of file Par03EMShowerModel.hh.

113{ return fNbOfHits; };

◆ SetLongMaxDepth()

void Par03EMShowerModel::SetLongMaxDepth ( const G4double  aDepth)
inline

Set maximum depth of shower created in fast simulation.

It is expressed in units of radiaton length.

Definition at line 116 of file Par03EMShowerModel.hh.

117 {
118 fLongMaxDepth = aDepth;
119 };

◆ GetLongMaxDepth()

G4double Par03EMShowerModel::GetLongMaxDepth ( ) const
inline

Get maximum depth of shower created in fast simulation.

It is expressed in units of radiaton length.

Definition at line 122 of file Par03EMShowerModel.hh.

122{ return fLongMaxDepth; };

◆ Gamma()

G4double Par03EMShowerModel::Gamma ( G4double  x,
G4double  alpha,
G4double  beta 
)
inlineprivate

Gamma distribution.

Definition at line 126 of file Par03EMShowerModel.hh.

127 {
128 return (std::pow(beta, alpha) / std::tgamma(alpha) * std::pow(x, alpha - 1) *
129 std::exp(-beta * x));
130 }

◆ Gaussian()

G4double Par03EMShowerModel::Gaussian ( G4double  x,
G4double  sigma = 1,
G4double  x0 = 0 
)
inlineprivate

Gaussian distribution.

Definition at line 132 of file Par03EMShowerModel.hh.

133 {
134 G4double tmp = (x - x0) / sigma;
135 return (1.0 / (std::sqrt(2 * CLHEP::pi) * sigma)) * std::exp(-tmp * tmp / 2);
136 }

Member Data Documentation

◆ fMessenger

Par03EMShowerMessenger* Par03EMShowerModel::fMessenger
private

Messenger for configuration.

Definition at line 140 of file Par03EMShowerModel.hh.

◆ fHitMaker

std::unique_ptr<G4FastSimHitMaker> Par03EMShowerModel::fHitMaker
private

Helper class for creation of hits within the sensitive detector.

Definition at line 142 of file Par03EMShowerModel.hh.

◆ fSigma

G4double Par03EMShowerModel::fSigma = -1
private

Standard deviation of the Gaussian distribution Can be changed with UI command /Par03/fastSim/transverseProfile/sigma <sigma> If sigma is smaller than 0, it will be estimated from the detector material (Moliere radius).

Definition at line 148 of file Par03EMShowerModel.hh.

◆ fAlpha

G4double Par03EMShowerModel::fAlpha = -1
private

Alpha parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudunalProfile/alpha <alpha> If alpha is smaller than 0, it will be estimated from particle energy and the detector material.

Definition at line 154 of file Par03EMShowerModel.hh.

◆ fBeta

G4double Par03EMShowerModel::fBeta = 0.5
private

Beta parameter of the Gamma distribution Can be changed with UI command /Par03/fastSim/longitudinalProfile/beta <beta>

Definition at line 158 of file Par03EMShowerModel.hh.

◆ fNbOfHits

G4int Par03EMShowerModel::fNbOfHits = 100
private

Number of (same energy) hits created by the parametrisation.

Can be changed with UI command /Par03/fastSim/numberOfHits <number>

Definition at line 161 of file Par03EMShowerModel.hh.

◆ fLongMaxDepth

G4double Par03EMShowerModel::fLongMaxDepth = 30
private

Maximum depth of a shower created in fast simulation.

It is expressed in units of radiation length. Can be changed with UI command /Par03/fastSim/longitudinalProfile/maxDepth <depth>

Definition at line 165 of file Par03EMShowerModel.hh.


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

Applications | User Support | Publications | Collaboration