Geant4 11.2.2
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+
67
68// hadrons
70#include "G4MscStepLimitType.hh"
71
72#include "G4ePairProduction.hh"
73
74#include "G4hIonisation.hh"
75#include "G4ionIonisation.hh"
77#include "G4NuclearStopping.hh"
79
80// msc models
81#include "G4UrbanMscModel.hh"
82#include "G4WentzelVIModel.hh"
86
87// interfaces
88#include "G4EmParameters.hh"
89#include "G4EmBuilder.hh"
90#include "G4EmStandUtil.hh"
91
92// particles
93
94#include "G4Gamma.hh"
95#include "G4Electron.hh"
96#include "G4Positron.hh"
97#include "G4GenericIon.hh"
98
99//
100#include "G4PhysicsListHelper.hh"
101#include "G4BuilderType.hh"
102#include "G4EmModelActivator.hh"
103
104// factory
106//
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112 : G4VPhysicsConstructor(pname)
113{
114 SetVerboseLevel(ver);
116 param->SetDefaults();
117 param->SetVerbose(ver);
118 param->SetMinEnergy(100*CLHEP::eV);
119 param->SetLowestElectronEnergy(100*CLHEP::eV);
120 param->SetNumberOfBinsPerDecade(20);
122 param->SetStepFunction(0.2, 10*CLHEP::um);
123 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
124 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
125 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
126 param->SetUseMottCorrection(true);
128 param->SetMscSkin(3);
129 param->SetMscRangeFactor(0.08);
130 param->SetMuHadLateralDisplacement(true);
131 param->SetFluo(true);
132 param->SetUseICRU90Data(true);
134 param->SetMaxNIELEnergy(1*CLHEP::MeV);
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
146{
147 // minimal set of particles for EM physics
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
154{
155 if(verboseLevel > 1) {
156 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
157 }
161
162 // processes used by several particles
163 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
164
165 // high energy limit for e+- scattering models
166 G4double highEnergyLimit= param->MscEnergyLimit();
167 G4double livEnergyLimit = 1*CLHEP::GeV;
168
169 // nuclear stopping
170 G4double nielEnergyLimit = param->MaxNIELEnergy();
171 G4NuclearStopping* pnuc = nullptr;
172 if(nielEnergyLimit > 0.0) {
173 pnuc = new G4NuclearStopping();
174 pnuc->SetMaxKinEnergy(nielEnergyLimit);
175 }
176
177 // Add Livermore EM Processes
178
179 // Add gamma EM Processes
181 G4bool polar = param->EnablePolarisation();
182
183 // photoelectric effect - Livermore model only
186 pe->SetEmModel(peModel);
187 if(polar) {
189 }
190
191 // Compton scattering - Livermore model only
194 G4VEmModel* cModel = nullptr;
195 if(polar) {
197 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
198 } else {
199 cModel = new G4LivermoreComptonModel();
200 cModel->SetHighEnergyLimit(livEnergyLimit);
201 }
202 cs->AddEmModel(0, cModel);
203
204 // gamma conversion
207 gc->SetEmModel(convLiv);
208
209 // Rayleigh scattering
211 if(polar) {
213 }
214
215 ph->RegisterProcess(pe, particle);
216 ph->RegisterProcess(cs, particle);
217 ph->RegisterProcess(gc, particle);
218 ph->RegisterProcess(rl, particle);
219
220 // e-
221 particle = G4Electron::Electron();
222
223 // multiple and single scattering
226 msc1->SetHighEnergyLimit(highEnergyLimit);
227 msc2->SetLowEnergyLimit(highEnergyLimit);
228 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
229
232 ss->SetEmModel(ssm);
233 ss->SetMinKinEnergy(highEnergyLimit);
234 ssm->SetLowEnergyLimit(highEnergyLimit);
235 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
236
237 // Ionisation - Livermore should be used only for low energies
238 G4eIonisation* eioni = new G4eIonisation();
240 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
241 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
242 eioni->AddEmModel(0, theIoniLiv);
243
244 // Bremsstrahlung
250 brem->SetEmModel(br1);
251 brem->SetEmModel(br2);
252 br1->SetHighEnergyLimit(GeV);
253
255
256 // register processes
257 ph->RegisterProcess(eioni, particle);
258 ph->RegisterProcess(brem, particle);
259 ph->RegisterProcess(ee, particle);
260 ph->RegisterProcess(ss, particle);
261
262 // e+
263 particle = G4Positron::Positron();
264
265 // multiple and single scattering
266 msc1 = new G4GoudsmitSaundersonMscModel();
267 msc2 = new G4WentzelVIModel();
268 msc1->SetHighEnergyLimit(highEnergyLimit);
269 msc2->SetLowEnergyLimit(highEnergyLimit);
270 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
271
272 ssm = new G4eCoulombScatteringModel();
273 ss = new G4CoulombScattering();
274 ss->SetEmModel(ssm);
275 ss->SetMinKinEnergy(highEnergyLimit);
276 ssm->SetLowEnergyLimit(highEnergyLimit);
277 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
278
279 // ionisation
280 eioni = new G4eIonisation();
281
282 // Bremsstrahlung from standard
283 brem = new G4eBremsstrahlung();
284 br1 = new G4SeltzerBergerModel();
285 br2 = new G4eBremsstrahlungRelModel();
288 brem->SetEmModel(br1);
289 brem->SetEmModel(br2);
290 br1->SetHighEnergyLimit(GeV);
291
292 // register processes
293 ph->RegisterProcess(eioni, particle);
294 ph->RegisterProcess(brem, particle);
295 ph->RegisterProcess(ee, particle);
296 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
297 ph->RegisterProcess(ss, particle);
298
299 // generic ion
300 particle = G4GenericIon::GenericIon();
301 G4ionIonisation* ionIoni = new G4ionIonisation();
303 ph->RegisterProcess(hmsc, particle);
304 ph->RegisterProcess(ionIoni, particle);
305 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
306
307 // muons, hadrons, ions
309
310 // extra configuration
312
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
@ 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()
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
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)
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)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)