Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmStandardPhysics.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
30//
31// Author: V.Ivanchenko 09.11.2005
32//
33// Modified:
34//
35//----------------------------------------------------------------------------
36//
37
39#include "G4SystemOfUnits.hh"
41#include "G4EmParameters.hh"
42#include "G4EmBuilder.hh"
43#include "G4LossTableManager.hh"
44
47#include "G4GammaConversion.hh"
53
57#include "G4WentzelVIModel.hh"
58#include "G4UrbanMscModel.hh"
59
60#include "G4eIonisation.hh"
61#include "G4eBremsstrahlung.hh"
63
64#include "G4hIonisation.hh"
65#include "G4ionIonisation.hh"
66#include "G4NuclearStopping.hh"
67
68#include "G4Gamma.hh"
69#include "G4Electron.hh"
70#include "G4Positron.hh"
71#include "G4GenericIon.hh"
72
74#include "G4BuilderType.hh"
75#include "G4EmModelActivator.hh"
77
78// factory
80//
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86 : G4VPhysicsConstructor("G4EmStandard")
87{
88 SetVerboseLevel(ver);
90 param->SetDefaults();
91 param->SetVerbose(ver);
92 param->SetGeneralProcessActive(true);
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
100{}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105{
106 // minimal set of particles for EM physics
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111
113{
114 if(verboseLevel > 1) {
115 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
116 }
118
121
122 // processes used by several particles
123 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
124
125 // nuclear stopping is enabled if th eenergy limit above zero
126 G4double nielEnergyLimit = param->MaxNIELEnergy();
127 G4NuclearStopping* pnuc = nullptr;
128 if(nielEnergyLimit > 0.0) {
129 pnuc = new G4NuclearStopping();
130 pnuc->SetMaxKinEnergy(nielEnergyLimit);
131 }
132
133 // high energy limit for e+- scattering models and bremsstrahlung
134 G4double highEnergyLimit = param->MscEnergyLimit();
135
136 // Add gamma EM Processes
138 G4bool polar = param->EnablePolarisation();
139
140 // Photoelectric
143 pe->SetEmModel(peModel);
144 if(polar) {
146 }
147
148 // Compton scattering
150 if(polar) {
152 }
153
154 // default Rayleigh scattering is Livermore
156 if(polar) {
158 }
159
160 if(G4EmParameters::Instance()->GeneralProcessActive()) {
162 sp->AddEmProcess(pe);
163 sp->AddEmProcess(cs);
164 sp->AddEmProcess(new G4GammaConversion());
165 sp->AddEmProcess(rl);
167 ph->RegisterProcess(sp, particle);
168
169 } else {
170 ph->RegisterProcess(pe, particle);
171 ph->RegisterProcess(cs, particle);
172 ph->RegisterProcess(new G4GammaConversion(), particle);
173 ph->RegisterProcess(rl, particle);
174 }
175
176 // e-
177 particle = G4Electron::Electron();
178
179 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
181 msc1->SetHighEnergyLimit(highEnergyLimit);
182 msc2->SetLowEnergyLimit(highEnergyLimit);
183 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
184
187 ss->SetEmModel(ssm);
188 ss->SetMinKinEnergy(highEnergyLimit);
189 ssm->SetLowEnergyLimit(highEnergyLimit);
190 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
191
192 ph->RegisterProcess(new G4eIonisation(), particle);
193 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
194 ph->RegisterProcess(ss, particle);
195
196 // e+
197 particle = G4Positron::Positron();
198
199 msc1 = new G4UrbanMscModel();
200 msc2 = new G4WentzelVIModel();
201 msc1->SetHighEnergyLimit(highEnergyLimit);
202 msc2->SetLowEnergyLimit(highEnergyLimit);
203 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
204
205 ssm = new G4eCoulombScatteringModel();
206 ss = new G4CoulombScattering();
207 ss->SetEmModel(ssm);
208 ss->SetMinKinEnergy(highEnergyLimit);
209 ssm->SetLowEnergyLimit(highEnergyLimit);
210 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
211
212 ph->RegisterProcess(new G4eIonisation(), particle);
213 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
214 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
215 ph->RegisterProcess(ss, particle);
216
217 // generic ion
218 particle = G4GenericIon::GenericIon();
219 G4ionIonisation* ionIoni = new G4ionIonisation();
220 ph->RegisterProcess(hmsc, particle);
221 ph->RegisterProcess(ionIoni, particle);
222 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
223
224 // muons, hadrons ions
226
227 // extra configuration
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
#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: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 ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
Definition: G4EmBuilder.cc:409
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:399
static G4EmParameters * Instance()
void SetGeneralProcessActive(G4bool val)
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
void SetFluctuationType(G4EmFluctuationType val)
void SetVerbose(G4int val)
G4EmStandardPhysics(G4int ver=1, const G4String &name="")
void ConstructParticle() override
void ConstructProcess() override
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: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)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)