Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmPenelopePhysics.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
29#include "G4SystemOfUnits.hh"
30
31// *** Processes and models
32
33// gamma
36
39
40#include "G4GammaConversion.hh"
42
45
48
49// e- and e+
51
52#include "G4eIonisation.hh"
54
55#include "G4eBremsstrahlung.hh"
57
58// e+ only
61
62// hadrons
64#include "G4MscStepLimitType.hh"
65
66#include "G4ePairProduction.hh"
67
68#include "G4hIonisation.hh"
69#include "G4ionIonisation.hh"
71#include "G4NuclearStopping.hh"
73
74// msc models
75#include "G4UrbanMscModel.hh"
77#include "G4WentzelVIModel.hh"
80
81// utilities
82#include "G4EmBuilder.hh"
83#include "G4EmStandUtil.hh"
84
85// particles
86#include "G4Gamma.hh"
87#include "G4Electron.hh"
88#include "G4Positron.hh"
89#include "G4GenericIon.hh"
90
91//
93#include "G4BuilderType.hh"
94#include "G4EmModelActivator.hh"
95
96// factory
98//
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104 : G4VPhysicsConstructor("G4EmPenelope")
105{
106 SetVerboseLevel(ver);
108 param->SetDefaults();
109 param->SetVerbose(ver);
110 param->SetMinEnergy(100*CLHEP::eV);
111 param->SetLowestElectronEnergy(100*CLHEP::eV);
112 param->SetNumberOfBinsPerDecade(20);
114 param->SetStepFunction(0.2, 10*CLHEP::um);
115 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
116 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
117 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
118 param->SetUseMottCorrection(true);
120 param->SetMscSkin(3);
121 param->SetMscRangeFactor(0.08);
122 param->SetMuHadLateralDisplacement(true);
123 param->SetFluo(true);
124 param->SetUseICRU90Data(true);
126 param->SetMaxNIELEnergy(1*CLHEP::MeV);
127 param->SetPIXEElectronCrossSectionModel("Penelope");
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
134{}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139{
140 // minimal set of particles for EM physics
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
147{
148 if(verboseLevel > 1) {
149 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
150 }
154
155 // processes used by several particles
156 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
157
158 // high energy limit for e+- scattering models
159 G4double highEnergyLimit = param->MscEnergyLimit();
160
161 // nuclear stopping
162 G4double nielEnergyLimit = param->MaxNIELEnergy();
163 G4NuclearStopping* pnuc = nullptr;
164 if(nielEnergyLimit > 0.0) {
165 pnuc = new G4NuclearStopping();
166 pnuc->SetMaxKinEnergy(nielEnergyLimit);
167 }
168
169 //Applicability range for Penelope models
170 //for higher energies, the Standard models are used
171 G4double PenelopeHighEnergyLimit = 1.0*CLHEP::GeV;
172
173 // Add gamma EM Processes
175
176 //Photo-electric effect
177 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
178 thePhotoElectricEffect->SetEmModel(new G4PEEffectFluoModel());
180 thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
181 thePhotoElectricEffect->AddEmModel(0, thePEPenelopeModel);
182
183 //Compton scattering
184 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
185 G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
186 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
187 theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
188 theComptonScattering->AddEmModel(0, theComptonPenelopeModel);
189
190 //Gamma conversion
191 G4GammaConversion* theGammaConversion = new G4GammaConversion();
193 theGCPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
194 theGammaConversion->AddEmModel(0, theGCPenelopeModel);
195
196 //Rayleigh scattering
197 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
198 G4PenelopeRayleighModel* theRayleighPenelopeModel = new G4PenelopeRayleighModel();
199 theRayleigh->SetEmModel(theRayleighPenelopeModel);
200
201 ph->RegisterProcess(thePhotoElectricEffect, particle);
202 ph->RegisterProcess(theComptonScattering, particle);
203 ph->RegisterProcess(theGammaConversion, particle);
204 ph->RegisterProcess(theRayleigh, particle);
205
206 // e-
207 particle = G4Electron::Electron();
208
209 // multiple scattering
212 msc1->SetHighEnergyLimit(highEnergyLimit);
213 msc2->SetLowEnergyLimit(highEnergyLimit);
214 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
215
218 ss->SetEmModel(ssm);
219 ss->SetMinKinEnergy(highEnergyLimit);
220 ssm->SetLowEnergyLimit(highEnergyLimit);
221 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
222
223 //Ionisation
224 G4eIonisation* eioni = new G4eIonisation();
227 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
228 eioni->AddEmModel(0, theIoniPenelope);
229
230 //Bremsstrahlung
233 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
234 brem->SetEmModel(theBremPenelope);
235
237
238 // register processes
239 ph->RegisterProcess(eioni, particle);
240 ph->RegisterProcess(brem, particle);
241 ph->RegisterProcess(ee, particle);
242 ph->RegisterProcess(ss, particle);
243
244 // e+
245 particle = G4Positron::Positron();
246
247 // multiple scattering
248 msc1 = new G4GoudsmitSaundersonMscModel();
249 msc2 = new G4WentzelVIModel();
250 msc1->SetHighEnergyLimit(highEnergyLimit);
251 msc2->SetLowEnergyLimit(highEnergyLimit);
252 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
253
254 ssm = new G4eCoulombScatteringModel();
255 ss = new G4CoulombScattering();
256 ss->SetEmModel(ssm);
257 ss->SetMinKinEnergy(highEnergyLimit);
258 ssm->SetLowEnergyLimit(highEnergyLimit);
259 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
260
261 //Ionisation
262 eioni = new G4eIonisation();
264 theIoniPenelope = new G4PenelopeIonisationModel();
265 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
266 eioni->AddEmModel(0,theIoniPenelope);
267
268 //Bremsstrahlung
269 brem = new G4eBremsstrahlung();
270 theBremPenelope = new G4PenelopeBremsstrahlungModel();
271 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
272 brem->SetEmModel(theBremPenelope);
273
274 //Annihilation
277 theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
278 anni->AddEmModel(0, theAnnPenelope);
279
280 // register processes
281 ph->RegisterProcess(eioni, particle);
282 ph->RegisterProcess(brem, particle);
283 ph->RegisterProcess(ee, particle);
284 ph->RegisterProcess(anni, particle);
285 ph->RegisterProcess(ss, particle);
286
287 // generic ion
288 particle = G4GenericIon::GenericIon();
289 G4ionIonisation* ionIoni = new G4ionIonisation();
291 ph->RegisterProcess(hmsc, particle);
292 ph->RegisterProcess(ionIoni, particle);
293 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
294
295 // muons, hadrons, ions
297
298 // extra configuration
300}
301
302//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseSafetyPlus
#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 ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
Definition: G4EmBuilder.cc:409
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()
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
void SetPIXEElectronCrossSectionModel(const G4String &)
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 SetMscSkin(G4double 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)
void ConstructParticle() override
G4EmPenelopePhysics(G4int ver=1, const G4String &name="")
void ConstructProcess() override
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
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 AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)