Geant4-11
G4EmLivermorePhysics.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// gamma
36
42
43#include "G4GammaConversion.hh"
46
50
53
54// e+-
57
58#include "G4eIonisation.hh"
60
61#include "G4eBremsstrahlung.hh"
63#include "G4Generator2BS.hh"
66
67// e+
69
70// hadrons
72#include "G4MscStepLimitType.hh"
73
74#include "G4ePairProduction.hh"
75
76#include "G4hIonisation.hh"
77#include "G4ionIonisation.hh"
79#include "G4NuclearStopping.hh"
80
81// msc models
82#include "G4UrbanMscModel.hh"
83#include "G4WentzelVIModel.hh"
87
88// interfaces
89#include "G4EmParameters.hh"
90#include "G4EmBuilder.hh"
91
92// particles
93
94#include "G4Gamma.hh"
95#include "G4Electron.hh"
96#include "G4Positron.hh"
97#include "G4GenericIon.hh"
98
99//
100#include "G4PhysicsListHelper.hh"
101#include "G4BuilderType.hh"
102#include "G4EmModelActivator.hh"
103
104// factory
106//
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
113{
114 SetVerboseLevel(ver);
116 param->SetDefaults();
117 param->SetVerbose(ver);
118 param->SetMinEnergy(100*CLHEP::eV);
120 param->SetNumberOfBinsPerDecade(20);
122 param->SetStepFunction(0.2, 10*CLHEP::um);
123 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
124 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
125 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
126 param->SetUseMottCorrection(true);
128 param->SetMscSkin(3);
129 param->SetMscRangeFactor(0.08);
130 param->SetMuHadLateralDisplacement(true);
131 param->SetFluo(true);
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139{}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
144{
145 // minimal set of particles for EM physics
147}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
152{
153 if(verboseLevel > 1) {
154 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
155 }
159
160 // processes used by several particles
161 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
162
163 // high energy limit for e+- scattering models
164 G4double highEnergyLimit= param->MscEnergyLimit();
165 G4double livEnergyLimit = 1*CLHEP::GeV;
166
167 // nuclear stopping
168 G4double nielEnergyLimit = param->MaxNIELEnergy();
169 G4NuclearStopping* pnuc = nullptr;
170 if(nielEnergyLimit > 0.0) {
171 pnuc = new G4NuclearStopping();
172 pnuc->SetMaxKinEnergy(nielEnergyLimit);
173 }
174
175 // Add Livermore EM Processes
176
177 // Add gamma EM Processes
179 G4bool polar = param->EnablePolarisation();
180
181 // photoelectric effect - Livermore model only
184 pe->SetEmModel(peModel);
185 if(polar) {
187 }
188
189 // Compton scattering - Livermore model only
192 G4VEmModel* cModel = nullptr;
193 if(polar) {
195 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
196 } else {
197 cModel = new G4LivermoreComptonModel();
198 cModel->SetHighEnergyLimit(livEnergyLimit);
199 }
200 cs->AddEmModel(0, cModel);
201
202 // gamma conversion
205 gc->SetEmModel(convLiv);
206
207 // Rayleigh scattering
209 if(polar) {
211 }
212
213 ph->RegisterProcess(pe, particle);
214 ph->RegisterProcess(cs, particle);
215 ph->RegisterProcess(gc, particle);
216 ph->RegisterProcess(rl, particle);
217
218 // e-
219 particle = G4Electron::Electron();
220
221 // multiple and single scattering
225 msc1->SetHighEnergyLimit(highEnergyLimit);
226 msc2->SetLowEnergyLimit(highEnergyLimit);
227 msc->SetEmModel(msc1);
228 msc->SetEmModel(msc2);
229
232 ss->SetEmModel(ssm);
233 ss->SetMinKinEnergy(highEnergyLimit);
234 ssm->SetLowEnergyLimit(highEnergyLimit);
235 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
236
237 // Ionisation - Livermore should be used only for low energies
238 G4eIonisation* eioni = new G4eIonisation();
239 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
240 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
241 eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation() );
242
243 // Bremsstrahlung
249 brem->SetEmModel(br1);
250 brem->SetEmModel(br2);
252
254
255 // register processes
256 ph->RegisterProcess(msc, particle);
257 ph->RegisterProcess(eioni, particle);
258 ph->RegisterProcess(brem, particle);
259 ph->RegisterProcess(ee, particle);
260 ph->RegisterProcess(ss, particle);
261
262 // e+
263 particle = G4Positron::Positron();
264
265 // multiple and single scattering
266 msc = new G4eMultipleScattering;
267 msc1 = new G4GoudsmitSaundersonMscModel();
268 msc2 = new G4WentzelVIModel();
269 msc1->SetHighEnergyLimit(highEnergyLimit);
270 msc2->SetLowEnergyLimit(highEnergyLimit);
271 msc->SetEmModel(msc1);
272 msc->SetEmModel(msc2);
273
274 ssm = new G4eCoulombScatteringModel();
275 ss = new G4CoulombScattering();
276 ss->SetEmModel(ssm);
277 ss->SetMinKinEnergy(highEnergyLimit);
278 ssm->SetLowEnergyLimit(highEnergyLimit);
279 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
280
281 // ionisation
282 eioni = new G4eIonisation();
283
284 // Bremsstrahlung from standard
285 brem = new G4eBremsstrahlung();
286 br1 = new G4SeltzerBergerModel();
287 br2 = new G4eBremsstrahlungRelModel();
290 brem->SetEmModel(br1);
291 brem->SetEmModel(br2);
293
294 // register processes
295 ph->RegisterProcess(msc, particle);
296 ph->RegisterProcess(eioni, particle);
297 ph->RegisterProcess(brem, particle);
298 ph->RegisterProcess(ee, particle);
299 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
300 ph->RegisterProcess(ss, particle);
301
302 // generic ion
303 particle = G4GenericIon::GenericIon();
304 G4ionIonisation* ionIoni = new G4ionIonisation();
306 ph->RegisterProcess(hmsc, particle);
307 ph->RegisterProcess(ionIoni, particle);
308 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
309
310 // muons, hadrons, ions
312
313 // extra configuration
315
316}
317
318//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLivermorePhysics)
@ fUseSafetyPlus
static constexpr double GeV
Definition: G4SIunits.hh:203
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:219
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:375
void ConstructProcess() override
void ConstructParticle() override
G4EmLivermorePhysics(G4int ver=1, const G4String &name="G4EmLivermore")
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() 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 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)
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 SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
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
string pname
Definition: eplot.py:33