Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmStandardPhysicsGS.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// ClassName: G4EmStandardPhysicsGS
30//
31// Author: V.Ivanchenko 05.06.2015
32//
33// Modified:
34//
35// 25.10.25 Changed the settings to the most accurate configuration of the
36// Goudsmit-Saunderson MSC model for e-/e+ multiple Coulomb
37// scattering. The class descrition has been changed to reflect the
38// new configuration (Mihaly Novak).
39//
40// Class Description:
41//
42// This EM physics constructor utilises the Goudsmit-Saunderson (GS) MSC model
43// for e-/e+ multiple Coulomb scattering (below 1 GeV kinetic energies). The
44// GS MSC model has been changed in version 11.4 keeping only its accurate
45// stepping and boundary crossing algorithms while removeing the other, less
46// accurate alternatives. All the corrections offered by the GS model, including
47// the Mott, screening and scattering power corrections, are activated. This,
48// together with the Penelope model for e-/e+ ionisations (below 1 GeV),
49// offers an accurate an accurate e-/e+ simualtion down to few keV kinetic
50// enegies independently form the target material and geometry.
51//
52// The same settings of the GS MSC model has already been used for e-/e+ below
53// 100 MeV kinetic energies in the option4, Penelope and Livermore EM physics
54// constructors since version 10.6.
55//
56//----------------------------------------------------------------------------
57//
58
60
61#include "G4SystemOfUnits.hh"
63#include "G4LossTableManager.hh"
64#include "G4EmParameters.hh"
65#include "G4EmStandUtil.hh"
66#include "G4EmBuilder.hh"
67
69#include "G4GammaConversion.hh"
73
77#include "G4WentzelVIModel.hh"
79
80#include "G4eIonisation.hh"
82#include "G4eBremsstrahlung.hh"
84#include "G4Generator2BS.hh"
85
87
88#include "G4hIonisation.hh"
89#include "G4ionIonisation.hh"
90
91#include "G4Gamma.hh"
92#include "G4Electron.hh"
93#include "G4Positron.hh"
94#include "G4GenericIon.hh"
95
97#include "G4BuilderType.hh"
98#include "G4EmModelActivator.hh"
100
101// factory
103//
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109 : G4VPhysicsConstructor("G4EmStandardGS")
110{
111 SetVerboseLevel(ver);
113 param->SetDefaults();
114 param->SetVerbose(ver);
115 // use a denser discrete kinetic energy grid for more accurate interpolation
116 param->SetNumberOfBinsPerDecade(16);
117 // set the continuous step limit to: 0.2*Range that goes to Range below 10 um
118 param->SetStepFunction(0.2, 10*CLHEP::um);
119 // set the GS MSC model for e-/e+ to be used below 1.0 GeV with its (Mott,
120 // screening, scattering power) corrections activated and with the accurate
121 // stepping and boundary crossing algorithms (no other options since 11.4)
122 // with a skin of 3 elastic MFP near boundary
123 param->SetMscEnergyLimit(1.0*CLHEP::GeV);
124 param->SetUseMottCorrection(true);
126 param->SetMscSkin(3);
127 param->SetMscRangeFactor(0.08);
128 // activate fluoresence, i.e. emission of characteristic X-ray
129 param->SetFluo(true);
130 // set the energy loss fluctuation type
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
143{
144 // minimal set of particles for EM physics
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
151{
152 if(verboseLevel > 1) {
153 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
154 }
157
158 // processes used by several particles
159 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
160 G4NuclearStopping* pnuc(nullptr);
161
162 // high energy limit for e+- scattering models and bremsstrahlung
164
165 // Add gamma EM processes
166
167 // gamma
169
172
176
177 if (G4EmParameters::Instance()->GeneralProcessActive()) {
179 sp->AddEmProcess(pe);
180 sp->AddEmProcess(cs);
181 sp->AddEmProcess(gc);
182 sp->AddEmProcess(rs);
184 ph->RegisterProcess(sp, particle);
185 } else {
186 ph->RegisterProcess(pe, particle);
187 ph->RegisterProcess(cs, particle);
188 ph->RegisterProcess(gc, particle);
189 ph->RegisterProcess(rs, particle);
190 }
191
192 // e-
193 particle = G4Electron::Electron();
194
195 // msc: GS[:100 MeV] + WentzelVI[100 MeV:]
198 msc1->SetHighEnergyLimit(mscEnergyLimit);
199 msc2->SetLowEnergyLimit(mscEnergyLimit);
200 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
201 // (WVI is a mixed model, i.e. needs single scattering)
204 ss->SetEmModel(ssm);
205 ss->SetMinKinEnergy(mscEnergyLimit);
206 ssm->SetLowEnergyLimit(mscEnergyLimit);
207 ssm->SetActivationLowEnergyLimit(mscEnergyLimit);
208
209 // ionisation: Penelope[:1.0 GeV] + Moller[1.0 GeV:]
210 G4eIonisation* eioni = new G4eIonisation();
212 G4VEmModel* theIoniMod = new G4PenelopeIonisationModel();
213 theIoniMod->SetHighEnergyLimit(1.0*CLHEP::GeV);
214 eioni->AddEmModel(0, theIoniMod);
215
216 // bremsstrahlung: Seltzer-Berger[:1.0 GeV] + extended Bethe–Heitler[1.0 GeV:]
222 brem->SetEmModel(br1);
223 brem->SetEmModel(br2);
224 br1->SetHighEnergyLimit(1.0*CLHEP::GeV);
225
226 ph->RegisterProcess(eioni, particle);
227 ph->RegisterProcess(brem, particle);
228 ph->RegisterProcess(ss, particle);
229
230 // e+
231 particle = G4Positron::Positron();
232
233 // msc: GS[:100 MeV] + WentzelVI[100 MeV:]
234 msc1 = new G4GoudsmitSaundersonMscModel();
235 msc2 = new G4WentzelVIModel();
236 msc1->SetHighEnergyLimit(mscEnergyLimit);
237 msc2->SetLowEnergyLimit(mscEnergyLimit);
238 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
239 // (WVI is a mixed model, i.e. needs single scattering)
240 ssm = new G4eCoulombScatteringModel();
241 ss = new G4CoulombScattering();
242 ss->SetEmModel(ssm);
243 ss->SetMinKinEnergy(mscEnergyLimit);
244 ssm->SetLowEnergyLimit(mscEnergyLimit);
245 ssm->SetActivationLowEnergyLimit(mscEnergyLimit);
246
247 // ionisation: Penelope[:1.0 GeV] + Bhabha[1.0 GeV:]
248 eioni = new G4eIonisation();
251 pen->SetHighEnergyLimit(1.0*CLHEP::GeV);
252 eioni->AddEmModel(0, pen);
253
254 // bremsstrahlung: Seltzer-Berger[:1.0 GeV] + extended Bethe–Heitler[1.0 GeV:]
255 brem = new G4eBremsstrahlung();
256 br1 = new G4SeltzerBergerModel();
257 br2 = new G4eBremsstrahlungRelModel();
260 brem->SetEmModel(br1);
261 brem->SetEmModel(br2);
262 br1->SetHighEnergyLimit(1.0*CLHEP::GeV);
263
264 ph->RegisterProcess(eioni, particle);
265 ph->RegisterProcess(brem, particle);
266 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
267 ph->RegisterProcess(ss, particle);
268
269 // generic ion
270 particle = G4GenericIon::GenericIon();
271 G4ionIonisation* ionIoni = new G4ionIonisation();
272 ph->RegisterProcess(hmsc, particle);
273 ph->RegisterProcess(ionIoni, particle);
274
275 // muons, hadrons, ions
277
278 // extra configuration
280}
281
282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseSafetyPlus
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
static void ConstructMinimalEmSet()
static void ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
G4double MscEnergyLimit() const
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetFluctuationType(G4EmFluctuationType val)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMscEnergyLimit(G4double val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
G4EmStandardPhysicsGS(G4int ver=1, const G4String &name="")
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
G4VPhysicsConstructor(const G4String &="")
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)