Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
79#include "G4IonFluctuations.hh"
80#include "G4NuclearStopping.hh"
81
82#include "G4ParticleTable.hh"
83#include "G4Gamma.hh"
84#include "G4Electron.hh"
85#include "G4Positron.hh"
86#include "G4GenericIon.hh"
87
89#include "G4BuilderType.hh"
90#include "G4EmModelActivator.hh"
92
93// factory
95//
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99
101 const G4String&)
102 : G4VPhysicsConstructor("G4EmStandard_opt3")
103{
104 SetVerboseLevel(ver);
106 param->SetDefaults();
107 param->SetVerbose(ver);
108 param->SetGeneralProcessActive(true);
109 param->SetMinEnergy(10*CLHEP::eV);
110 param->SetLowestElectronEnergy(100*CLHEP::eV);
111 param->SetNumberOfBinsPerDecade(20);
113 param->SetUseMottCorrection(true);
114 param->SetStepFunction(0.2, 100*CLHEP::um);
115 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
116 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
117 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
119 param->SetMuHadLateralDisplacement(true);
120 param->SetLateralDisplacementAlg96(true);
121 param->SetUseICRU90Data(true);
123 param->SetFluo(true);
124 param->SetMaxNIELEnergy(1*CLHEP::MeV);
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136{
137 // minimal set of particles for EM physics
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
144{
145 if(verboseLevel > 1) {
146 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
147 }
149
152
153 // processes used by several particles
154 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
155
156 // nuclear stopping is enabled if th eenergy limit above zero
157 G4double nielEnergyLimit = param->MaxNIELEnergy();
158 G4NuclearStopping* pnuc = nullptr;
159 if(nielEnergyLimit > 0.0) {
160 pnuc = new G4NuclearStopping();
161 pnuc->SetMaxKinEnergy(nielEnergyLimit);
162 }
163
164 // Add gamma EM Processes
166
169 pe->SetEmModel(peModel);
170 if(param->EnablePolarisation()) {
172 }
173
176
178 if(param->EnablePolarisation()) {
180 }
181
183 if(param->EnablePolarisation()) {
185 }
186
187 if(G4EmParameters::Instance()->GeneralProcessActive()) {
189 sp->AddEmProcess(pe);
190 sp->AddEmProcess(cs);
191 sp->AddEmProcess(gc);
192 sp->AddEmProcess(rl);
194 ph->RegisterProcess(sp, particle);
195 } else {
196 ph->RegisterProcess(pe, particle);
197 ph->RegisterProcess(cs, particle);
198 ph->RegisterProcess(gc, particle);
199 ph->RegisterProcess(rl, particle);
200 }
201
202 // e-
203 particle = G4Electron::Electron();
204
205 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
206 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
207
208 G4eIonisation* eIoni = new G4eIonisation();
209
215 brem->SetEmModel(br1);
216 brem->SetEmModel(br2);
217 br2->SetLowEnergyLimit(CLHEP::GeV);
218
220
221 ph->RegisterProcess(eIoni, particle);
222 ph->RegisterProcess(brem, particle);
223 ph->RegisterProcess(ee, particle);
224
225 // e+
226 particle = G4Positron::Positron();
227
228 msc1 = new G4UrbanMscModel();
229 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
230
231 eIoni = new G4eIonisation();
232
233 brem = new G4eBremsstrahlung();
234 br1 = new G4SeltzerBergerModel();
235 br2 = new G4eBremsstrahlungRelModel();
238 brem->SetEmModel(br1);
239 brem->SetEmModel(br2);
240 br2->SetLowEnergyLimit(CLHEP::GeV);
241
242 ph->RegisterProcess(eIoni, particle);
243 ph->RegisterProcess(brem, particle);
244 ph->RegisterProcess(ee, particle);
245 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
246
247 // generic ion
248 particle = G4GenericIon::GenericIon();
249 G4ionIonisation* ionIoni = new G4ionIonisation();
250 auto fluc = new G4IonFluctuations();
251 ionIoni->SetFluctModel(fluc);
253 ph->RegisterProcess(hmsc, particle);
254 ph->RegisterProcess(ionIoni, particle);
255 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
256
257 // muons, hadrons, ions
258 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
259
260 // extra configuration
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseDistanceToBoundary
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition G4Types.hh:83
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
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 SetFluctuationType(G4EmFluctuationType 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)
G4EmStandardPhysics_option3(G4int ver=1, const G4String &name="")
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 SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)