Geant4-11
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
34
37
40
41#include "G4GammaConversion.hh"
43
46
49
50// e- and e+
51
54
55#include "G4eIonisation.hh"
57
58#include "G4eBremsstrahlung.hh"
60
61// e+ only
62
65
66// hadrons
67
69#include "G4MscStepLimitType.hh"
70
71#include "G4ePairProduction.hh"
72
73#include "G4hIonisation.hh"
74#include "G4ionIonisation.hh"
76#include "G4NuclearStopping.hh"
77
78// msc models
79#include "G4UrbanMscModel.hh"
81#include "G4WentzelVIModel.hh"
84//
85
86#include "G4EmBuilder.hh"
87
88// particles
89
90#include "G4Gamma.hh"
91#include "G4Electron.hh"
92#include "G4Positron.hh"
93#include "G4GenericIon.hh"
94
95//
97#include "G4BuilderType.hh"
98#include "G4EmModelActivator.hh"
99
100// factory
102//
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
108 : G4VPhysicsConstructor("G4EmPenelope")
109{
110 SetVerboseLevel(ver);
112 param->SetDefaults();
113 param->SetVerbose(ver);
114 param->SetMinEnergy(100*CLHEP::eV);
116 param->SetNumberOfBinsPerDecade(20);
118 param->SetStepFunction(0.2, 10*CLHEP::um);
119 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
120 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
121 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
122 param->SetUseMottCorrection(true);
124 param->SetMscSkin(3);
125 param->SetMscRangeFactor(0.08);
126 param->SetMuHadLateralDisplacement(true);
127 param->SetFluo(true);
129 param->SetPIXEElectronCrossSectionModel("Penelope");
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136{}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141{
142 // minimal set of particles for EM physics
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
149{
150 if(verboseLevel > 1) {
151 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
152 }
156
157 // processes used by several particles
158 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
159
160 // high energy limit for e+- scattering models
161 G4double highEnergyLimit = param->MscEnergyLimit();
162
163 // nuclear stopping
164 G4double nielEnergyLimit = param->MaxNIELEnergy();
165 G4NuclearStopping* pnuc = nullptr;
166 if(nielEnergyLimit > 0.0) {
167 pnuc = new G4NuclearStopping();
168 pnuc->SetMaxKinEnergy(nielEnergyLimit);
169 }
170
171 //Applicability range for Penelope models
172 //for higher energies, the Standard models are used
173 G4double PenelopeHighEnergyLimit = 1.0*CLHEP::GeV;
174
175 // Add gamma EM Processes
177
178 //Photo-electric effect
179 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
180 thePhotoElectricEffect->SetEmModel(new G4PEEffectFluoModel());
182 thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
183 thePhotoElectricEffect->AddEmModel(0, thePEPenelopeModel);
184
185 //Compton scattering
186 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
187 G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
188 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
189 theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
190 theComptonScattering->AddEmModel(0, theComptonPenelopeModel);
191
192 //Gamma conversion
193 G4GammaConversion* theGammaConversion = new G4GammaConversion();
195 theGCPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
196 theGammaConversion->AddEmModel(0, theGCPenelopeModel);
197
198 //Rayleigh scattering
199 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
200 G4PenelopeRayleighModel* theRayleighPenelopeModel = new G4PenelopeRayleighModel();
201 theRayleigh->SetEmModel(theRayleighPenelopeModel);
202
203 ph->RegisterProcess(thePhotoElectricEffect, particle);
204 ph->RegisterProcess(theComptonScattering, particle);
205 ph->RegisterProcess(theGammaConversion, particle);
206 ph->RegisterProcess(theRayleigh, particle);
207
208 // e-
209 particle = G4Electron::Electron();
210
211 // multiple scattering
215 msc1->SetHighEnergyLimit(highEnergyLimit);
216 msc2->SetLowEnergyLimit(highEnergyLimit);
217 msc->SetEmModel(msc1);
218 msc->SetEmModel(msc2);
219
222 ss->SetEmModel(ssm);
223 ss->SetMinKinEnergy(highEnergyLimit);
224 ssm->SetLowEnergyLimit(highEnergyLimit);
225 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
226
227 //Ionisation
228 G4eIonisation* eioni = new G4eIonisation();
230 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
231 eioni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
232
233 //Bremsstrahlung
236 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
237 brem->SetEmModel(theBremPenelope);
238
240
241 // register processes
242 ph->RegisterProcess(msc, particle);
243 ph->RegisterProcess(eioni, particle);
244 ph->RegisterProcess(brem, particle);
245 ph->RegisterProcess(ee, particle);
246 ph->RegisterProcess(ss, particle);
247
248 // e+
249 particle = G4Positron::Positron();
250
251 // multiple scattering
252 msc = new G4eMultipleScattering;
253 msc1 = new G4GoudsmitSaundersonMscModel();
254 msc2 = new G4WentzelVIModel();
255 msc1->SetHighEnergyLimit(highEnergyLimit);
256 msc2->SetLowEnergyLimit(highEnergyLimit);
257 msc->SetEmModel(msc1);
258 msc->SetEmModel(msc2);
259
260 ssm = new G4eCoulombScatteringModel();
261 ss = new G4CoulombScattering();
262 ss->SetEmModel(ssm);
263 ss->SetMinKinEnergy(highEnergyLimit);
264 ssm->SetLowEnergyLimit(highEnergyLimit);
265 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
266
267 //Ionisation
268 eioni = new G4eIonisation();
269 theIoniPenelope = new G4PenelopeIonisationModel();
270 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
271 eioni->AddEmModel(0,theIoniPenelope, new G4UniversalFluctuation());
272
273 //Bremsstrahlung
274 brem = new G4eBremsstrahlung();
275 theBremPenelope = new G4PenelopeBremsstrahlungModel();
276 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
277 brem->SetEmModel(theBremPenelope);
278
279 //Annihilation
282 theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
283 anni->AddEmModel(0, theAnnPenelope);
284
285 // register processes
286 ph->RegisterProcess(msc, particle);
287 ph->RegisterProcess(eioni, particle);
288 ph->RegisterProcess(brem, particle);
289 ph->RegisterProcess(ee, particle);
290 ph->RegisterProcess(anni, particle);
291 ph->RegisterProcess(ss, particle);
292
293 // generic ion
294 particle = G4GenericIon::GenericIon();
295 G4ionIonisation* ionIoni = new G4ionIonisation();
297 ph->RegisterProcess(hmsc, particle);
298 ph->RegisterProcess(ionIoni, particle);
299 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
300
301 // muons, hadrons, ions
303
304 // extra configuration
306}
307
308//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmPenelopePhysics)
@ fUseSafetyPlus
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:219
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:375
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 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 SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
void ConstructParticle() override
G4EmPenelopePhysics(G4int ver=1, const G4String &name="")
void ConstructProcess() override
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:767
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
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 SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)
static constexpr double um
Definition: SystemOfUnits.h:94
static constexpr double GeV
static constexpr double MeV
static constexpr double eV