Geant4 11.1.1
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// | History: |
27// | -------- |
28// | |
29// | Feb. 2021 JMCB - Adapted for polarized gamma ray transport. |
30// | See "An electromagnetic physics constructor |
31// | for low energy polarised X-/gamma ray transport |
32// | in Geant4", J. M. C. Brown and M. R. Dimmock, |
33// | arXiv:2102.02721 (2021). |
34// | https://arxiv.org/abs/2102.02721 |
35// | |
36// *********************************************************************
37//
38
39#include "G4EmLowEPPhysics.hh"
41#include "G4SystemOfUnits.hh"
42#include "G4ParticleTable.hh"
43
44// *** Processes and models
45
46// gamma
50
53
54#include "G4GammaConversion.hh"
56
59
62
63// e+-
66#include "G4ePairProduction.hh"
67
68#include "G4eIonisation.hh"
70
71#include "G4eBremsstrahlung.hh"
73#include "G4Generator2BS.hh"
74
75// e+
78
79// hadrons
81#include "G4MscStepLimitType.hh"
82
83#include "G4hIonisation.hh"
84#include "G4ionIonisation.hh"
87#include "G4NuclearStopping.hh"
88
89// msc models
90#include "G4UrbanMscModel.hh"
91#include "G4WentzelVIModel.hh"
96
97// interfaces
98#include "G4LossTableManager.hh"
99#include "G4EmBuilder.hh"
100#include "G4EmParameters.hh"
101
102// particles
103
104#include "G4Gamma.hh"
105#include "G4Electron.hh"
106#include "G4Positron.hh"
107#include "G4GenericIon.hh"
108
109//
110#include "G4PhysicsListHelper.hh"
111#include "G4BuilderType.hh"
112#include "G4EmModelActivator.hh"
113
114// factory
116//
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
122 : G4VPhysicsConstructor("G4EmLowEPPhysics")
123{
124 SetVerboseLevel(ver);
126 param->SetDefaults();
127 param->SetVerbose(ver);
128 param->SetMinEnergy(100*CLHEP::eV);
129 param->SetLowestElectronEnergy(100*CLHEP::eV);
130 param->SetNumberOfBinsPerDecade(20);
132 param->SetStepFunction(0.2, 100*CLHEP::um);
133 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
134 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
135 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
136 param->SetUseMottCorrection(true);
137 param->SetMscRangeFactor(0.04);
138 param->SetMuHadLateralDisplacement(true);
139 param->SetFluo(true);
140 param->SetAuger(true);
141 param->SetUseICRU90Data(true);
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
148{}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
153{
154 // minimal set of particles for EM physics
156}
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159
161{
162 if(verboseLevel > 1) {
163 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
164 }
166
169
170 // common processes
171 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
172
173 // nuclear stopping
174 G4double nielEnergyLimit = param->MaxNIELEnergy();
175 G4NuclearStopping* pnuc = nullptr;
176 if(nielEnergyLimit > 0.0) {
177 pnuc = new G4NuclearStopping();
178 pnuc->SetMaxKinEnergy(nielEnergyLimit);
179 }
180
181 // high energy limit for e+- scattering models and bremsstrahlung
182 G4double highEnergyLimit = param->MscEnergyLimit();
183
184 // Add gamma EM Processes
186
187 // Photoelectric absorption
189 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
191 pe->SetEmModel(theLivermorePEModel);
192
193 // Compton scattering - Polarised Monash model
195 G4VEmModel* theLowEPCSModel = new G4LowEPPolarizedComptonModel();
196 cs->SetEmModel(theLowEPCSModel);
197
198 // gamma conversion - 5D model below 80 GeV with Livermore x-sections
199 G4GammaConversion* theGammaConversion = new G4GammaConversion();
200 G4VEmModel* conv = new G4BetheHeitler5DModel();
201 theGammaConversion->SetEmModel(conv);
202
203 // default Rayleigh scattering is Livermore
204 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
205 G4VEmModel* theLivermorePRSModel = new G4LivermorePolarizedRayleighModel();
206 theRayleigh->SetEmModel(theLivermorePRSModel);
207
208 ph->RegisterProcess(pe, particle);
209 ph->RegisterProcess(cs, particle);
210 ph->RegisterProcess(theGammaConversion, particle);
211 ph->RegisterProcess(theRayleigh, particle);
212
213 // e-
214 particle = G4Electron::Electron();
215
216 // multiple scattering
220 msc1->SetHighEnergyLimit(highEnergyLimit);
221 msc2->SetLowEnergyLimit(highEnergyLimit);
222 msc->SetEmModel(msc1);
223 msc->SetEmModel(msc2);
224
227 ss->SetEmModel(ssm);
228 ss->SetMinKinEnergy(highEnergyLimit);
229 ssm->SetLowEnergyLimit(highEnergyLimit);
230 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
231
232 // Ionisation - Livermore should be used only for low energies
233 G4eIonisation* eioni = new G4eIonisation();
235 theIoniLivermore->SetHighEnergyLimit(0.1*CLHEP::MeV);
236 eioni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
237
238 // Bremsstrahlung
244 brem->SetEmModel(br1);
245 brem->SetEmModel(br2);
246 br1->SetHighEnergyLimit(GeV);
247
249
250 // register processes
251 ph->RegisterProcess(msc, particle);
252 ph->RegisterProcess(ss, particle);
253 ph->RegisterProcess(eioni, particle);
254 ph->RegisterProcess(brem, particle);
255 ph->RegisterProcess(ee, particle);
256
257 // e+
258 particle = G4Positron::Positron();
259
260 // multiple scattering
261 msc = new G4eMultipleScattering();
262 msc1 = new G4GoudsmitSaundersonMscModel();
263 msc2 = new G4WentzelVIModel();
264 msc1->SetHighEnergyLimit(highEnergyLimit);
265 msc2->SetLowEnergyLimit(highEnergyLimit);
266 msc->SetEmModel(msc1);
267 msc->SetEmModel(msc2);
268
269 ssm = new G4eCoulombScatteringModel();
270 ss = new G4CoulombScattering();
271 ss->SetEmModel(ssm);
272 ss->SetMinKinEnergy(highEnergyLimit);
273 ssm->SetLowEnergyLimit(highEnergyLimit);
274 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
275
276 // Standard ionisation
277 eioni = new G4eIonisation();
279 pen->SetHighEnergyLimit(0.1*MeV);
280 eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
281
282 // Bremsstrahlung
283 brem = new G4eBremsstrahlung();
284 br1 = new G4SeltzerBergerModel();
285 br2 = new G4eBremsstrahlungRelModel();
288 brem->SetEmModel(br1);
289 brem->SetEmModel(br2);
290 br1->SetHighEnergyLimit(CLHEP::GeV);
291
292 // register processes
293 ph->RegisterProcess(msc, particle);
294 ph->RegisterProcess(ss, particle);
295 ph->RegisterProcess(eioni, particle);
296 ph->RegisterProcess(brem, particle);
297 ph->RegisterProcess(ee, particle);
298 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
299
300 // generic ion
301 particle = G4GenericIon::GenericIon();
302
303 G4ionIonisation* ionIoni = new G4ionIonisation();
305 ionIoni->SetEmModel(mod);
306
307 ph->RegisterProcess(hmsc, particle);
308 ph->RegisterProcess(ionIoni, particle);
309 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
310
311 // muons, hadrons ions
313
314 // extra configuration
316}
317
318//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition: G4Types.hh:83
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:232
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:399
~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()
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 SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetAuger(G4bool val)
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:746
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)