Geant4 10.7.0
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"
55
58#include "G4MscStepLimitType.hh"
59#include "G4UrbanMscModel.hh"
60#include "G4DummyModel.hh"
61#include "G4WentzelVIModel.hh"
63
64#include "G4eIonisation.hh"
65#include "G4eBremsstrahlung.hh"
66#include "G4Generator2BS.hh"
68
71
72#include "G4ePairProduction.hh"
73#include "G4ionIonisation.hh"
75#include "G4NuclearStopping.hh"
76
77#include "G4ParticleTable.hh"
78#include "G4Gamma.hh"
79#include "G4Electron.hh"
80#include "G4Positron.hh"
81#include "G4GenericIon.hh"
82
84#include "G4BuilderType.hh"
85#include "G4EmModelActivator.hh"
87
88// factory
90//
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94
96 const G4String&)
97 : G4VPhysicsConstructor("G4EmStandard_opt3"), verbose(ver)
98{
100 param->SetDefaults();
101 param->SetVerbose(verbose);
102 param->SetMinEnergy(10*CLHEP::eV);
103 param->SetLowestElectronEnergy(100*CLHEP::eV);
104 param->SetNumberOfBinsPerDecade(20);
106 param->SetUseMottCorrection(true);
107 param->SetStepFunction(0.2, 100*CLHEP::um);
108 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
109 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
110 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
112 param->SetMuHadLateralDisplacement(true);
113 param->SetLateralDisplacementAlg96(true);
114 param->SetUseICRU90Data(true);
115 param->SetFluo(true);
116 param->SetMaxNIELEnergy(1*CLHEP::MeV);
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121
123{}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
128{
129 // minimal set of particles for EM physics
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136{
137 if(verbose > 1) {
138 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
139 }
141
143
144 // processes used by several particles
146 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
147
148 // nuclear stopping
149 G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
151 pnuc->SetMaxKinEnergy(nielEnergyLimit);
152
153 // Add gamma EM Processes
155
158
161
162 if(G4EmParameters::Instance()->GeneralProcessActive()) {
164 sp->AddEmProcess(pee);
165 sp->AddEmProcess(cs);
166 sp->AddEmProcess(new G4GammaConversion());
167 sp->AddEmProcess(new G4RayleighScattering());
169 ph->RegisterProcess(sp, particle);
170 } else {
171 ph->RegisterProcess(pee,particle);
172 ph->RegisterProcess(cs, particle);
173 ph->RegisterProcess(new G4GammaConversion(), particle);
174 ph->RegisterProcess(new G4RayleighScattering(), particle);
175 }
176
177 // e-
178 particle = G4Electron::Electron();
179
181 G4eIonisation* eIoni = new G4eIonisation();
182
188 brem->SetEmModel(br1);
189 brem->SetEmModel(br2);
190 br2->SetLowEnergyLimit(CLHEP::GeV);
191
192 ph->RegisterProcess(msc, particle);
193 ph->RegisterProcess(eIoni, particle);
194 ph->RegisterProcess(brem, particle);
195 ph->RegisterProcess(ee, particle);
196
197 // e+
198 particle = G4Positron::Positron();
199
200 msc = new G4eMultipleScattering();
201 eIoni = new G4eIonisation();
202
203 brem = new G4eBremsstrahlung();
204 br1 = new G4SeltzerBergerModel();
205 br2 = new G4eBremsstrahlungRelModel();
208 brem->SetEmModel(br1);
209 brem->SetEmModel(br2);
210 br2->SetLowEnergyLimit(CLHEP::GeV);
211
212 ph->RegisterProcess(msc, particle);
213 ph->RegisterProcess(eIoni, particle);
214 ph->RegisterProcess(brem, particle);
215 ph->RegisterProcess(ee, particle);
216 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
217
218 // generic ion
219 particle = G4GenericIon::GenericIon();
220 G4ionIonisation* ionIoni = new G4ionIonisation();
222 ph->RegisterProcess(hmsc, particle);
223 ph->RegisterProcess(ionIoni, particle);
224 ph->RegisterProcess(pnuc, particle);
225
226 // muons, hadrons, ions
227 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
228
229 // extra configuration
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ 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: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:170
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
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)
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:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const