Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmLowEPPhysics.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#include "G4EmLowEPPhysics.hh"
29#include "G4SystemOfUnits.hh"
30#include "G4ParticleTable.hh"
31
32// *** Processes and models
33
34// gamma
37
41
42#include "G4GammaConversion.hh"
44
47
50
51// e+-
54#include "G4ePairProduction.hh"
55
56#include "G4eIonisation.hh"
58
59#include "G4eBremsstrahlung.hh"
61#include "G4Generator2BS.hh"
63
64// e+
66
67// hadrons
69#include "G4MscStepLimitType.hh"
70
71#include "G4hIonisation.hh"
72#include "G4ionIonisation.hh"
75#include "G4NuclearStopping.hh"
76
77// msc models
78#include "G4UrbanMscModel.hh"
79#include "G4WentzelVIModel.hh"
84
85// interfaces
86#include "G4LossTableManager.hh"
87#include "G4EmBuilder.hh"
88#include "G4EmParameters.hh"
89
90// particles
91
92#include "G4Gamma.hh"
93#include "G4Electron.hh"
94#include "G4Positron.hh"
95#include "G4GenericIon.hh"
96
97//
99#include "G4BuilderType.hh"
100#include "G4EmModelActivator.hh"
101
102// factory
104//
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
110 : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
111{
113 param->SetDefaults();
114 param->SetVerbose(verbose);
115 param->SetMinEnergy(100*CLHEP::eV);
116 param->SetLowestElectronEnergy(100*CLHEP::eV);
117 param->SetNumberOfBinsPerDecade(20);
119 param->SetStepFunction(0.2, 100*CLHEP::um);
120 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
121 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
122 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
123 param->SetUseMottCorrection(true);
124 param->SetMscRangeFactor(0.04);
125 param->SetMuHadLateralDisplacement(true);
126 param->SetFluo(true);
127 param->SetUseICRU90Data(true);
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
134{}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139{
140 // minimal set of particles for EM physics
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
147{
148 if(verbose > 1) {
149 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
150 }
152
154
155 // common processes
157 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
158
159 // nuclear stopping
161 pnuc->SetMaxKinEnergy(CLHEP::MeV);
162
163 // Add gamma EM Processes
165
167 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
168 pe->SetEmModel(theLivermorePEModel);
169
170 // Compton scattering - Livermore model above 20 MeV, Monarsh's model below
173 G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
174 theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
175 cs->AddEmModel(0, theLowEPComptonModel);
176
177 // gamma conversion - 5D model below 80 GeV with Livermore x-sections
178 G4GammaConversion* theGammaConversion = new G4GammaConversion();
179 G4VEmModel* conv = new G4BetheHeitler5DModel();
180 theGammaConversion->SetEmModel(conv);
181
182 // default Rayleigh scattering is Livermore
183 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
184
185 ph->RegisterProcess(pe, particle);
186 ph->RegisterProcess(cs, particle);
187 ph->RegisterProcess(theGammaConversion, particle);
188 ph->RegisterProcess(theRayleigh, particle);
189
190 // e-
191 particle = G4Electron::Electron();
192
193 // multiple scattering
196
197 // Ionisation - Livermore should be used only for low energies
198 G4eIonisation* eioni = new G4eIonisation();
200 theIoniLivermore->SetHighEnergyLimit(0.1*CLHEP::MeV);
201 eioni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
202
203 // Bremsstrahlung
209 brem->SetEmModel(br1);
210 brem->SetEmModel(br2);
211 br1->SetHighEnergyLimit(GeV);
212
213 // register processes
214 ph->RegisterProcess(msc, particle);
215 ph->RegisterProcess(eioni, particle);
216 ph->RegisterProcess(brem, particle);
217
218 // e+
219 particle = G4Positron::Positron();
220
221 // multiple scattering
222 msc = new G4eMultipleScattering();
224
225 // Standard ionisation
226 eioni = new G4eIonisation();
227
228 // Bremsstrahlung
229 brem = new G4eBremsstrahlung();
230 br1 = new G4SeltzerBergerModel();
231 br2 = new G4eBremsstrahlungRelModel();
234 brem->SetEmModel(br1);
235 brem->SetEmModel(br2);
236 br1->SetHighEnergyLimit(CLHEP::GeV);
237
238 // register processes
239 ph->RegisterProcess(msc, particle);
240 ph->RegisterProcess(eioni, particle);
241 ph->RegisterProcess(brem, particle);
242 ph->RegisterProcess(ee, particle);
243 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
244
245 // generic ion
246 particle = G4GenericIon::GenericIon();
247
248 G4ionIonisation* ionIoni = new G4ionIonisation();
250 ionIoni->SetEmModel(mod);
251
252 ph->RegisterProcess(hmsc, particle);
253 ph->RegisterProcess(ionIoni, particle);
254 ph->RegisterProcess(pnuc, particle);
255
256 // muons, hadrons ions
258
259 // extra configuration
261}
262
263//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
Definition: G4EmBuilder.cc:170
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
~G4EmLowEPPhysics() override
G4EmLowEPPhysics(G4int ver=1, const G4String &name="")
void ConstructProcess() override
void ConstructParticle() override
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetEmModel(G4VMscModel *, size_t index=0)
const G4String & GetPhysicsName() const