Geant4-11
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
51#include "G4GammaConversion.hh"
62
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"
82#include "G4UrbanFluctuation.hh"
83
85
86#include "G4hIonisation.hh"
87#include "G4ionIonisation.hh"
89#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->SetMinEnergy(100*CLHEP::eV);
118 param->SetNumberOfBinsPerDecade(20);
120 param->SetStepFunction(0.2, 10*CLHEP::um);
121 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
122 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
123 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
124 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
125 param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
126 param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
127 param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
128 param->SetMuHadLateralDisplacement(true);
129 param->SetFluo(true);
130 param->SetUseICRU90Data(true);
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136
138{}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
143{
144 // minimal set of particles for EM physics
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
151{
152 if(verboseLevel > 1) {
153 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
154 }
156
159
160 // processes used by several particles
161 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
162
163 // nuclear stopping is enabled if the energy limit above zero
164 G4double nielEnergyLimit = param->MaxNIELEnergy();
165 G4NuclearStopping* pnuc = nullptr;
166 if(nielEnergyLimit > 0.0) {
167 pnuc = new G4NuclearStopping();
168 pnuc->SetMaxKinEnergy(nielEnergyLimit);
169 }
170
171 // high energy limit for e+- scattering models and bremsstrahlung
172 G4double highEnergyLimit = param->MscEnergyLimit();
173
174 // Add gamma EM Processes
176 G4bool polar = param->EnablePolarisation();
177
178 // Photoelectric
181 pe->SetEmModel(peModel);
182 if(polar) {
184 }
185
186 // Compton scattering
189 G4VEmModel* cModel = nullptr;
190 if(polar) {
191 cModel = new G4LowEPPolarizedComptonModel();
192 } else {
193 cModel = new G4LowEPComptonModel();
194 }
195 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
196 cs->AddEmModel(0, cModel);
197
198 // Gamma conversion
200 G4VEmModel* conv = new G4BetheHeitler5DModel();
201 gc->SetEmModel(conv);
202
203 // default Rayleigh scattering is Livermore
205 if(polar) {
207 }
208
209 if(param->GeneralProcessActive()) {
211 sp->AddEmProcess(pe);
212 sp->AddEmProcess(cs);
213 sp->AddEmProcess(gc);
214 sp->AddEmProcess(rl);
216 ph->RegisterProcess(sp, particle);
217 } else {
218 ph->RegisterProcess(pe, particle);
219 ph->RegisterProcess(cs, particle);
220 ph->RegisterProcess(gc, particle);
221 ph->RegisterProcess(rl, particle);
222 }
223
224 // e-
225 particle = G4Electron::Electron();
226
227 // multiple scattering
229 // e-/e+ msc gs with Mott-correction
230 // (Mott-correction is set through G4EmParameters)
233 msc1->SetHighEnergyLimit(highEnergyLimit);
234 msc2->SetLowEnergyLimit(highEnergyLimit);
235 msc->SetEmModel(msc1);
236 msc->SetEmModel(msc2);
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();
247 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
248 eioni->SetFluctModel(new G4UrbanFluctuation());
249 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
250 eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation() );
251
252 // bremsstrahlung
258 brem->SetEmModel(br1);
259 brem->SetEmModel(br2);
261
263
264 // register processes
265 ph->RegisterProcess(msc, particle);
266 ph->RegisterProcess(eioni, particle);
267 ph->RegisterProcess(brem, particle);
268 ph->RegisterProcess(ee, particle);
269 ph->RegisterProcess(ss, particle);
270
271 // e+
272 particle = G4Positron::Positron();
273
274 // multiple scattering
275 msc = new G4eMultipleScattering();
276 // e-/e+ msc gs with Mott-correction
277 // (Mott-correction is set through G4EmParameters)
278 msc1 = new G4GoudsmitSaundersonMscModel();
279 msc2 = new G4WentzelVIModel();
280 msc1->SetHighEnergyLimit(highEnergyLimit);
281 msc2->SetLowEnergyLimit(highEnergyLimit);
282 msc->SetEmModel(msc1);
283 msc->SetEmModel(msc2);
284
285 ssm = new G4eCoulombScatteringModel();
286 ss = new G4CoulombScattering();
287 ss->SetEmModel(ssm);
288 ss->SetMinKinEnergy(highEnergyLimit);
289 ssm->SetLowEnergyLimit(highEnergyLimit);
290 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
291
292 // ionisation
293 eioni = new G4eIonisation();
294 eioni->SetFluctModel(new G4UrbanFluctuation());
297 eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
298
299 // bremsstrahlung
300 brem = new G4eBremsstrahlung();
301 br1 = new G4SeltzerBergerModel();
302 br2 = new G4eBremsstrahlungRelModel();
305 brem->SetEmModel(br1);
306 brem->SetEmModel(br2);
308
309 // register processes
310 ph->RegisterProcess(msc, particle);
311 ph->RegisterProcess(eioni, particle);
312 ph->RegisterProcess(brem, particle);
313 ph->RegisterProcess(ee, particle);
314 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
315 ph->RegisterProcess(ss, particle);
316
317 // generic ion
318 particle = G4GenericIon::GenericIon();
319 G4ionIonisation* ionIoni = new G4ionIonisation();
321 ph->RegisterProcess(hmsc, particle);
322 ph->RegisterProcess(ionIoni, particle);
323 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
324
325 // muons, hadrons, ions
327
328 // extra configuration
330}
331
332//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmStandardPhysics_option4)
@ fUseSafetyPlus
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
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:219
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:375
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 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)
G4bool GeneralProcessActive() const
G4EmStandardPhysics_option4(G4int ver=1, const G4String &name="")
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
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)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)
static constexpr double um
Definition: SystemOfUnits.h:94
static constexpr double GeV
static constexpr double MeV
static constexpr double eV