Geant4-11
G4EmStandardPhysics_option3.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_option3
30//
31// Author: V.Ivanchenko 13.03.2008
32//
33// Modified:
34// 21.04.2008 V.Ivanchenko add long-lived D and B mesons; use spline
35// 28.05.2008 V.Ivanchenko linLossLimit=0.01 for ions 0.001 for others
36//
37//----------------------------------------------------------------------------
38//
39
41
42#include "G4SystemOfUnits.hh"
44#include "G4LossTableManager.hh"
45#include "G4EmParameters.hh"
46#include "G4EmBuilder.hh"
47
49#include "G4GammaConversion.hh"
58
61#include "G4MscStepLimitType.hh"
62#include "G4UrbanMscModel.hh"
63#include "G4DummyModel.hh"
64#include "G4WentzelVIModel.hh"
66
67#include "G4eIonisation.hh"
68#include "G4eBremsstrahlung.hh"
69#include "G4Generator2BS.hh"
71
74
75#include "G4ePairProduction.hh"
76#include "G4ionIonisation.hh"
78#include "G4NuclearStopping.hh"
79
80#include "G4ParticleTable.hh"
81#include "G4Gamma.hh"
82#include "G4Electron.hh"
83#include "G4Positron.hh"
84#include "G4GenericIon.hh"
85
87#include "G4BuilderType.hh"
88#include "G4EmModelActivator.hh"
90
91// factory
93//
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99 const G4String&)
100 : G4VPhysicsConstructor("G4EmStandard_opt3")
101{
102 SetVerboseLevel(ver);
104 param->SetDefaults();
105 param->SetVerbose(ver);
106 param->SetMinEnergy(10*CLHEP::eV);
108 param->SetNumberOfBinsPerDecade(20);
110 param->SetUseMottCorrection(true);
111 param->SetStepFunction(0.2, 100*CLHEP::um);
112 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
113 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
114 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
116 param->SetMscRangeFactor(0.03);
117 param->SetMuHadLateralDisplacement(true);
118 param->SetLateralDisplacementAlg96(true);
119 param->SetUseICRU90Data(true);
120 param->SetFluo(true);
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
128{}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131
133{
134 // minimal set of particles for EM physics
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141{
142 if(verboseLevel > 1) {
143 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
144 }
146
149
150 // processes used by several particles
151 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
152
153 // nuclear stopping is enabled if th eenergy limit above zero
154 G4double nielEnergyLimit = param->MaxNIELEnergy();
155 G4NuclearStopping* pnuc = nullptr;
156 if(nielEnergyLimit > 0.0) {
157 pnuc = new G4NuclearStopping();
158 pnuc->SetMaxKinEnergy(nielEnergyLimit);
159 }
160
161 // Add gamma EM Processes
163
166 pe->SetEmModel(peModel);
167 if(param->EnablePolarisation()) {
169 }
170
173
175 if(param->EnablePolarisation()) {
177 }
178
180 if(param->EnablePolarisation()) {
182 }
183
184 if(G4EmParameters::Instance()->GeneralProcessActive()) {
186 sp->AddEmProcess(pe);
187 sp->AddEmProcess(cs);
188 sp->AddEmProcess(gc);
189 sp->AddEmProcess(rl);
191 ph->RegisterProcess(sp, particle);
192 } else {
193 ph->RegisterProcess(pe, particle);
194 ph->RegisterProcess(cs, particle);
195 ph->RegisterProcess(gc, particle);
196 ph->RegisterProcess(rl, particle);
197 }
198
199 // e-
200 particle = G4Electron::Electron();
201
203 G4eIonisation* eIoni = new G4eIonisation();
204
210 brem->SetEmModel(br1);
211 brem->SetEmModel(br2);
213
215
216 ph->RegisterProcess(msc, particle);
217 ph->RegisterProcess(eIoni, particle);
218 ph->RegisterProcess(brem, particle);
219 ph->RegisterProcess(ee, particle);
220
221 // e+
222 particle = G4Positron::Positron();
223
224 msc = new G4eMultipleScattering();
225 eIoni = new G4eIonisation();
226
227 brem = new G4eBremsstrahlung();
228 br1 = new G4SeltzerBergerModel();
229 br2 = new G4eBremsstrahlungRelModel();
232 brem->SetEmModel(br1);
233 brem->SetEmModel(br2);
235
236 ph->RegisterProcess(msc, particle);
237 ph->RegisterProcess(eIoni, particle);
238 ph->RegisterProcess(brem, particle);
239 ph->RegisterProcess(ee, particle);
240 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
241
242 // generic ion
243 particle = G4GenericIon::GenericIon();
244 G4ionIonisation* ionIoni = new G4ionIonisation();
246 ph->RegisterProcess(hmsc, particle);
247 ph->RegisterProcess(ionIoni, particle);
248 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
249
250 // muons, hadrons, ions
251 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
252
253 // extra configuration
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmStandardPhysics_option3)
@ fUseSafetyPlus
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: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
void SetLateralDisplacementAlg96(G4bool val)
G4double MaxNIELEnergy() 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 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)
G4EmStandardPhysics_option3(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 SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=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