Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmStandardPhysicsWVI.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: G4EmStandardPhysicsWVI
30//
31// Author: V.Ivanchenko 09.11.2005
32//
33// Modified:
34//
35// 21.04.2008 V.Ivanchenko add long-lived D and B mesons
36//
37//----------------------------------------------------------------------------
38//
39
41#include "G4SystemOfUnits.hh"
43#include "G4EmParameters.hh"
44#include "G4EmBuilder.hh"
45#include "G4LossTableManager.hh"
46
48#include "G4GammaConversion.hh"
51
57#include "G4WentzelVIModel.hh"
59#include "G4UrbanMscModel.hh"
62
63#include "G4eIonisation.hh"
64#include "G4eBremsstrahlung.hh"
67
68#include "G4hIonisation.hh"
69#include "G4ionIonisation.hh"
75#include "G4BraggIonModel.hh"
76#include "G4IonFluctuations.hh"
77#include "G4NuclearStopping.hh"
79
80#include "G4Gamma.hh"
81#include "G4Electron.hh"
82#include "G4Positron.hh"
83#include "G4GenericIon.hh"
84
86#include "G4BuilderType.hh"
87#include "G4EmModelActivator.hh"
88
89// factory
91//
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
97 : G4VPhysicsConstructor("G4EmStandardWVI")
98{
100 param->SetDefaults();
101 param->SetVerbose(ver);
102 param->SetMinEnergy(10*CLHEP::eV);
103 param->SetLowestElectronEnergy(100*CLHEP::eV);
104 param->SetNumberOfBinsPerDecade(20);
106 param->SetStepFunction(0.2, 100*CLHEP::um);
107 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
108 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
109 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
110 param->SetUseMottCorrection(true);
111 param->SetMuHadLateralDisplacement(true);
112 param->SetUseICRU90Data(true);
113 param->SetMscThetaLimit(0.15);
114 param->SetFluo(true);
115 param->SetMaxNIELEnergy(1*CLHEP::MeV);
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
122{}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127{
128 // minimal set of particles for EM physics
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
135{
136 if(verboseLevel > 1) {
137 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
138 }
142
143 // common processes
144 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
145
146 // nuclear stopping is enabled if th eenergy limit above zero
147 G4double nielEnergyLimit = param->MaxNIELEnergy();
148 G4NuclearStopping* pnuc = nullptr;
149 if(nielEnergyLimit > 0.0) {
150 pnuc = new G4NuclearStopping();
151 pnuc->SetMaxKinEnergy(nielEnergyLimit);
152 }
153
154 // high energy limit for e+- scattering models
155 G4double highEnergyLimit = 1*CLHEP::MeV;
156
157 // Add gamma EM processes
159
162
165
167 if(param->EnablePolarisation()) {
169 }
170
171 ph->RegisterProcess(pee, particle);
172 ph->RegisterProcess(cs, particle);
173 ph->RegisterProcess(gc, particle);
174 ph->RegisterProcess(new G4RayleighScattering(), particle);
175
176 // e-
177 particle = G4Electron::Electron();
178
180 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
182 msc1->SetHighEnergyLimit(highEnergyLimit);
183 msc2->SetLowEnergyLimit(highEnergyLimit);
184 msc->SetEmModel(msc1);
185 msc->SetEmModel(msc2);
186
189 ss->SetEmModel(ssm);
190 ss->SetMinKinEnergy(highEnergyLimit);
191 ssm->SetLowEnergyLimit(highEnergyLimit);
192 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
193
194 ph->RegisterProcess(msc, particle);
195 ph->RegisterProcess(new G4eIonisation(), particle);
196 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
197 ph->RegisterProcess(ss, particle);
198
199 // e+
200 particle = G4Positron::Positron();
201
202 msc = new G4eMultipleScattering;
203 msc1 = new G4UrbanMscModel();
204 msc2 = new G4WentzelVIModel();
205 msc1->SetHighEnergyLimit(highEnergyLimit);
206 msc2->SetLowEnergyLimit(highEnergyLimit);
207 msc->SetEmModel(msc1);
208 msc->SetEmModel(msc2);
209
210 ssm = new G4eCoulombScatteringModel();
211 ss = new G4CoulombScattering();
212 ss->SetEmModel(ssm);
213 ss->SetMinKinEnergy(highEnergyLimit);
214 ssm->SetLowEnergyLimit(highEnergyLimit);
215 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
216
219
220 ph->RegisterProcess(msc, particle);
221 ph->RegisterProcess(new G4eIonisation(), particle);
222 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
223 ph->RegisterProcess(ann, particle);
224 ph->RegisterProcess(ss, particle);
225
226 // generic ion
227 particle = G4GenericIon::GenericIon();
228 G4ionIonisation* ionIoni = new G4ionIonisation();
229 auto fluc = new G4IonFluctuations();
230 ionIoni->SetFluctModel(fluc);
232 ph->RegisterProcess(hmsc, particle);
233 ph->RegisterProcess(ionIoni, particle);
234 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
235
236 // muons, hadrons, ions
238
239 // extra configuration
241}
242
243//....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
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
void SetMscThetaLimit(G4double val)
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 SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool 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 SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const