Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmStandardPhysics_option4.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: G4EmStandardPhysics_option4
30//
31// Author: V.Ivanchenko 28.09.2012 from Option3 physics constructor
32//
33// Modified:
34//
35// 22.09.17 M.Novak: change msc model for e-/e+ below 100 MeV from Urban+
36// UseDistanceToBoundary stepping to GS + Mott-correction + error
37// -free stepping.
38//
39//----------------------------------------------------------------------------
40//
41
43
44#include "G4SystemOfUnits.hh"
46#include "G4LossTableManager.hh"
47#include "G4EmParameters.hh"
48#include "G4EmBuilder.hh"
49#include "G4EmStandUtil.hh"
50
52#include "G4GammaConversion.hh"
63
65#include "G4MscStepLimitType.hh"
66#include "G4UrbanMscModel.hh"
68#include "G4DummyModel.hh"
69#include "G4WentzelVIModel.hh"
72
73#include "G4eIonisation.hh"
74#include "G4eBremsstrahlung.hh"
75#include "G4Generator2BS.hh"
76#include "G4Generator2BN.hh"
78#include "G4ePairProduction.hh"
81
83
84#include "G4hIonisation.hh"
85#include "G4ionIonisation.hh"
88#include "G4NuclearStopping.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 const G4String&)
110 : G4VPhysicsConstructor("G4EmStandard_opt4")
111{
112 SetVerboseLevel(ver);
114 param->SetDefaults();
115 param->SetVerbose(ver);
116 param->SetGeneralProcessActive(true);
117 param->SetMinEnergy(100*CLHEP::eV);
118 param->SetLowestElectronEnergy(100*CLHEP::eV);
119 param->SetNumberOfBinsPerDecade(20);
121 param->SetStepFunction(0.2, 10*CLHEP::um);
122 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
123 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
124 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
125 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
126 param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
127 param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
128 param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
129 param->SetMuHadLateralDisplacement(true);
130 param->SetFluo(true);
131 param->SetUseICRU90Data(true);
132 param->Set3GammaAnnihilationOnFly(true);
134 param->SetMaxNIELEnergy(1*CLHEP::MeV);
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
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 }
159
162
163 // processes used by several particles
164 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
165
166 // nuclear stopping is enabled if the energy limit above zero
167 G4double nielEnergyLimit = param->MaxNIELEnergy();
168 G4NuclearStopping* pnuc = nullptr;
169 if(nielEnergyLimit > 0.0) {
170 pnuc = new G4NuclearStopping();
171 pnuc->SetMaxKinEnergy(nielEnergyLimit);
172 }
173
174 // high energy limit for e+- scattering models and bremsstrahlung
175 G4double highEnergyLimit = param->MscEnergyLimit();
176
177 // Add gamma EM Processes
179 G4bool polar = param->EnablePolarisation();
180
181 // Photoelectric
184 pe->SetEmModel(peModel);
185 if(polar) {
187 }
188
189 // Compton scattering
192 G4VEmModel* cModel = nullptr;
193 if(polar) {
194 cModel = new G4LowEPPolarizedComptonModel();
195 } else {
196 cModel = new G4LowEPComptonModel();
197 }
198 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
199 cs->AddEmModel(0, cModel);
200
201 // Gamma conversion
203 G4VEmModel* conv = new G4BetheHeitler5DModel();
204 gc->SetEmModel(conv);
205
206 // default Rayleigh scattering is Livermore
208 if(polar) {
210 }
211
212 if(param->GeneralProcessActive()) {
214 sp->AddEmProcess(pe);
215 sp->AddEmProcess(cs);
216 sp->AddEmProcess(gc);
217 sp->AddEmProcess(rl);
219 ph->RegisterProcess(sp, particle);
220 } else {
221 ph->RegisterProcess(pe, particle);
222 ph->RegisterProcess(cs, particle);
223 ph->RegisterProcess(gc, particle);
224 ph->RegisterProcess(rl, particle);
225 }
226
227 // e-
228 particle = G4Electron::Electron();
229
230 // e-/e+ msc gs with Mott-correction
231 // (Mott-correction is set through G4EmParameters)
234 msc1->SetHighEnergyLimit(highEnergyLimit);
235 msc2->SetLowEnergyLimit(highEnergyLimit);
236 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
237
240 ss->SetEmModel(ssm);
241 ss->SetMinKinEnergy(highEnergyLimit);
242 ssm->SetLowEnergyLimit(highEnergyLimit);
243 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
244
245 // ionisation
246 G4eIonisation* eioni = new G4eIonisation();
248 G4VEmModel* theIoniMod = new G4PenelopeIonisationModel();
249 theIoniMod->SetHighEnergyLimit(0.1*CLHEP::MeV);
250 eioni->AddEmModel(0, theIoniMod);
251
252 // bremsstrahlung
258 brem->SetEmModel(br1);
259 brem->SetEmModel(br2);
260 br1->SetHighEnergyLimit(CLHEP::GeV);
261
263
264 // register processes
265 ph->RegisterProcess(eioni, particle);
266 ph->RegisterProcess(brem, particle);
267 ph->RegisterProcess(ee, particle);
268 ph->RegisterProcess(ss, particle);
269
270 // e+
271 particle = G4Positron::Positron();
272
273 // e-/e+ msc gs with Mott-correction
274 // (Mott-correction is set through G4EmParameters)
275 msc1 = new G4GoudsmitSaundersonMscModel();
276 msc2 = new G4WentzelVIModel();
277 msc1->SetHighEnergyLimit(highEnergyLimit);
278 msc2->SetLowEnergyLimit(highEnergyLimit);
279 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
280
281 ssm = new G4eCoulombScatteringModel();
282 ss = new G4CoulombScattering();
283 ss->SetEmModel(ssm);
284 ss->SetMinKinEnergy(highEnergyLimit);
285 ssm->SetLowEnergyLimit(highEnergyLimit);
286 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
287
288 // ionisation
289 eioni = new G4eIonisation();
292 pen->SetHighEnergyLimit(0.1*CLHEP::MeV);
293 eioni->AddEmModel(0, pen);
294
295 // bremsstrahlung
296 brem = new G4eBremsstrahlung();
297 br1 = new G4SeltzerBergerModel();
298 br2 = new G4eBremsstrahlungRelModel();
301 brem->SetEmModel(br1);
302 brem->SetEmModel(br2);
303 br1->SetHighEnergyLimit(CLHEP::GeV);
304
305 // annihilation
306 auto anni = new G4eplusAnnihilation();
307 if (param->Use3GammaAnnihilationOnFly()) {
308 anni->SetEmModel(new G4eplusTo2or3GammaModel());
309 }
310
311 // register processes
312 ph->RegisterProcess(eioni, particle);
313 ph->RegisterProcess(brem, particle);
314 ph->RegisterProcess(ee, particle);
315 ph->RegisterProcess(anni, particle);
316 ph->RegisterProcess(ss, particle);
317
318 // generic ion
319 particle = G4GenericIon::GenericIon();
320 G4ionIonisation* ionIoni = new G4ionIonisation();
322 ph->RegisterProcess(hmsc, particle);
323 ph->RegisterProcess(ionIoni, particle);
324 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
325
326 // muons, hadrons, ions
328
329 // extra configuration
331}
332
333//....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()
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetGeneralProcessActive(G4bool val)
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)
void Set3GammaAnnihilationOnFly(G4bool v)
G4bool GeneralProcessActive() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
G4EmStandardPhysics_option4(G4int ver=1, const G4String &name="")
~G4EmStandardPhysics_option4() override
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 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)