Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmLivermorePhysics.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
29#include "G4SystemOfUnits.hh"
30
31// *** Processes and models
32// gamma
36
42
43#include "G4GammaConversion.hh"
46
50
53
54// e+-
56#include "G4eIonisation.hh"
58
59#include "G4eBremsstrahlung.hh"
61#include "G4Generator2BS.hh"
64
65// e+
68
69// hadrons
71#include "G4MscStepLimitType.hh"
72
73#include "G4ePairProduction.hh"
74
75#include "G4hIonisation.hh"
76#include "G4ionIonisation.hh"
78#include "G4NuclearStopping.hh"
80
81// msc models
82#include "G4UrbanMscModel.hh"
83#include "G4WentzelVIModel.hh"
87
88// interfaces
89#include "G4EmParameters.hh"
90#include "G4EmBuilder.hh"
91#include "G4EmStandUtil.hh"
92
93// particles
94
95#include "G4Gamma.hh"
96#include "G4Electron.hh"
97#include "G4Positron.hh"
98#include "G4GenericIon.hh"
99
100//
101#include "G4PhysicsListHelper.hh"
102#include "G4BuilderType.hh"
103#include "G4EmModelActivator.hh"
104
105// factory
107//
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111
113 : G4VPhysicsConstructor(pname)
114{
115 SetVerboseLevel(ver);
117 param->SetDefaults();
118 param->SetVerbose(ver);
119 param->SetMinEnergy(100*CLHEP::eV);
120 param->SetLowestElectronEnergy(100*CLHEP::eV);
121 param->SetNumberOfBinsPerDecade(20);
123 param->SetStepFunction(0.2, 10*CLHEP::um);
124 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
125 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
126 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
127 param->SetUseMottCorrection(true);
129 param->SetMscSkin(3);
130 param->SetMscRangeFactor(0.08);
131 param->SetMuHadLateralDisplacement(true);
132 param->SetFluo(true);
133 param->SetUseICRU90Data(true);
135 param->SetMaxNIELEnergy(1*CLHEP::MeV);
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
147{
148 // minimal set of particles for EM physics
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153
155{
156 if(verboseLevel > 1) {
157 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
158 }
162
163 // processes used by several particles
164 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
165
166 // high energy limit for e+- scattering models
167 G4double highEnergyLimit= param->MscEnergyLimit();
168 G4double livEnergyLimit = 1*CLHEP::GeV;
169
170 // nuclear stopping
171 G4double nielEnergyLimit = param->MaxNIELEnergy();
172 G4NuclearStopping* pnuc = nullptr;
173 if(nielEnergyLimit > 0.0) {
174 pnuc = new G4NuclearStopping();
175 pnuc->SetMaxKinEnergy(nielEnergyLimit);
176 }
177
178 // Add Livermore EM Processes
179
180 // Add gamma EM Processes
182 G4bool polar = param->EnablePolarisation();
183
184 // photoelectric effect - Livermore model only
187 pe->SetEmModel(peModel);
188 if(polar) {
190 }
191
192 // Compton scattering - Livermore model only
195 G4VEmModel* cModel = nullptr;
196 if(polar) {
198 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
199 } else {
200 cModel = new G4LivermoreComptonModel();
201 cModel->SetHighEnergyLimit(livEnergyLimit);
202 }
203 cs->AddEmModel(0, cModel);
204
205 // gamma conversion
208 gc->SetEmModel(convLiv);
209
210 // Rayleigh scattering
212 if(polar) {
214 }
215
216 ph->RegisterProcess(pe, particle);
217 ph->RegisterProcess(cs, particle);
218 ph->RegisterProcess(gc, particle);
219 ph->RegisterProcess(rl, particle);
220
221 // e-
222 particle = G4Electron::Electron();
223
224 // multiple and single scattering
227 msc1->SetHighEnergyLimit(highEnergyLimit);
228 msc2->SetLowEnergyLimit(highEnergyLimit);
229 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
230
233 ss->SetEmModel(ssm);
234 ss->SetMinKinEnergy(highEnergyLimit);
235 ssm->SetLowEnergyLimit(highEnergyLimit);
236 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
237
238 // Ionisation - Livermore should be used only for low energies
239 G4eIonisation* eioni = new G4eIonisation();
241 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
242 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
243 eioni->AddEmModel(0, theIoniLiv);
244
245 // Bremsstrahlung
251 brem->SetEmModel(br1);
252 brem->SetEmModel(br2);
253 br1->SetHighEnergyLimit(GeV);
254
256
257 // register processes
258 ph->RegisterProcess(eioni, particle);
259 ph->RegisterProcess(brem, particle);
260 ph->RegisterProcess(ee, particle);
261 ph->RegisterProcess(ss, particle);
262
263 // e+
264 particle = G4Positron::Positron();
265
266 // multiple and single scattering
267 msc1 = new G4GoudsmitSaundersonMscModel();
268 msc2 = new G4WentzelVIModel();
269 msc1->SetHighEnergyLimit(highEnergyLimit);
270 msc2->SetLowEnergyLimit(highEnergyLimit);
271 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
272
273 ssm = new G4eCoulombScatteringModel();
274 ss = new G4CoulombScattering();
275 ss->SetEmModel(ssm);
276 ss->SetMinKinEnergy(highEnergyLimit);
277 ssm->SetLowEnergyLimit(highEnergyLimit);
278 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
279
280 // ionisation
281 eioni = new G4eIonisation();
282
283 // Bremsstrahlung from standard
284 brem = new G4eBremsstrahlung();
285 br1 = new G4SeltzerBergerModel();
286 br2 = new G4eBremsstrahlungRelModel();
289 brem->SetEmModel(br1);
290 brem->SetEmModel(br2);
291 br1->SetHighEnergyLimit(GeV);
292
293 // annihilation
294 auto anni = new G4eplusAnnihilation();
295 if (param->Use3GammaAnnihilationOnFly()) {
296 anni->SetEmModel(new G4eplusTo2or3GammaModel());
297 }
298
299 // register processes
300 ph->RegisterProcess(eioni, particle);
301 ph->RegisterProcess(brem, particle);
302 ph->RegisterProcess(ee, particle);
303 ph->RegisterProcess(anni, particle);
304 ph->RegisterProcess(ss, particle);
305
306 // generic ion
307 particle = G4GenericIon::GenericIon();
308 G4ionIonisation* ionIoni = new G4ionIonisation();
310 ph->RegisterProcess(hmsc, particle);
311 ph->RegisterProcess(ionIoni, particle);
312 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
313
314 // muons, hadrons, ions
316
317 // extra configuration
319
320}
321
322//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
@ fAllisonPositronium
@ fUseSafetyPlus
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
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()
~G4EmLivermorePhysics() override
void ConstructParticle() override
G4EmLivermorePhysics(G4int ver=1, const G4String &name="G4EmLivermore")
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
void SetPositronAtRestModelType(G4PositronAtRestModelType val)
G4double MscEnergyLimit() const
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetFluctuationType(G4EmFluctuationType val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
G4bool Use3GammaAnnihilationOnFly() const
void SetMscRangeFactor(G4double val)
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
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 AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
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)