Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmStandardPhysics_option2.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_option2
30//
31// Author: V.Ivanchenko 09.11.2005
32//
33// Modified:
34// 19.12.2005 V.Ivanchenko rename 71 -> 72
35// 15.06.2006 V.Ivanchenko use this class as a constructor of fast EM physics
36// 13.11.2006 V.Ivanchenko use G4hMultipleScattering
37// 14.11.2006 V.Ivanchenko use sub-cutoff option for all particles
38// 13.02.2007 V.Ivanchenko use default msc
39// 15.05.2007 V.Ivanchenko rename to _option2
40// 13.03.2008 V.Ivanchenko use G4eMultipleScattering
41// 21.04.2008 V.Ivanchenko add long-lived D and B mesons; use spline
42// 28.05.2008 V.Ivanchenko linLossLimit=0.01; added hBrem and hPairProd processes
43//
44//----------------------------------------------------------------------------
45//
46
48
49#include "G4SystemOfUnits.hh"
51#include "G4LossTableManager.hh"
52#include "G4EmParameters.hh"
53#include "G4EmBuilder.hh"
54
56#include "G4GammaConversion.hh"
61
66#include "G4UrbanMscModel.hh"
67#include "G4WentzelVIModel.hh"
68
69#include "G4eIonisation.hh"
70#include "G4eBremsstrahlung.hh"
73#include "G4Generator2BS.hh"
75
76#include "G4hIonisation.hh"
77#include "G4ionIonisation.hh"
78
79#include "G4ParticleTable.hh"
80#include "G4Gamma.hh"
81#include "G4Electron.hh"
82#include "G4Positron.hh"
83#include "G4GenericIon.hh"
84
86#include "G4BuilderType.hh"
87#include "G4EmModelActivator.hh"
89
90// factory
92//
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96
98 const G4String&)
99 : G4VPhysicsConstructor("G4EmStandard_opt2"), verbose(ver)
100{
102 param->SetDefaults();
103 param->SetVerbose(verbose);
104 param->SetApplyCuts(false);
105 param->SetStepFunction(0.8, 1*CLHEP::mm);
106 param->SetMscRangeFactor(0.2);
107 param->SetLateralDisplacement(false);
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113
115{}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
120{
121 // minimal set of particles for EM physics
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
128{
129 if(verbose > 1) {
130 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
131 }
133
135
136 // processes used by several particles
137 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
138 G4NuclearStopping* pnuc(nullptr);
139
140 // high energy limit for e+- scattering models and bremsstrahlung
141 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
142
143 // Add gamma EM Processes
145
148
149 if(G4EmParameters::Instance()->GeneralProcessActive()) {
151 sp->AddEmProcess(pee);
152 sp->AddEmProcess(new G4ComptonScattering());
153 sp->AddEmProcess(new G4GammaConversion());
155 ph->RegisterProcess(sp, particle);
156
157 } else {
158 ph->RegisterProcess(pee, particle);
159 ph->RegisterProcess(new G4ComptonScattering(), particle);
160 ph->RegisterProcess(new G4GammaConversion(), particle);
161 }
162
163 // e-
164 particle = G4Electron::Electron();
165
166 G4eIonisation* eioni = new G4eIonisation();
167
169 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
171 msc1->SetHighEnergyLimit(highEnergyLimit);
172 msc2->SetLowEnergyLimit(highEnergyLimit);
173 msc->SetEmModel(msc1);
174 msc->SetEmModel(msc2);
175
178 ss->SetEmModel(ssm);
179 ss->SetMinKinEnergy(highEnergyLimit);
180 ssm->SetLowEnergyLimit(highEnergyLimit);
181 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
182
183 ph->RegisterProcess(msc, particle);
184 ph->RegisterProcess(eioni, particle);
185 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
186 ph->RegisterProcess(ss, particle);
187
188 // e+
189 particle = G4Positron::Positron();
190 eioni = new G4eIonisation();
191
192 msc = new G4eMultipleScattering;
193 msc1 = new G4UrbanMscModel();
194 msc2 = new G4WentzelVIModel();
195 msc1->SetHighEnergyLimit(highEnergyLimit);
196 msc2->SetLowEnergyLimit(highEnergyLimit);
197 msc->SetEmModel(msc1);
198 msc->SetEmModel(msc2);
199
200 ssm = new G4eCoulombScatteringModel();
201 ss = new G4CoulombScattering();
202 ss->SetEmModel(ssm);
203 ss->SetMinKinEnergy(highEnergyLimit);
204 ssm->SetLowEnergyLimit(highEnergyLimit);
205 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
206
207 ph->RegisterProcess(msc, particle);
208 ph->RegisterProcess(eioni, particle);
209 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
210 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
211 ph->RegisterProcess(ss, particle);
212
213 // generic ion
214 particle = G4GenericIon::GenericIon();
215 G4ionIonisation* ionIoni = new G4ionIonisation();
216 ph->RegisterProcess(hmsc, particle);
217 ph->RegisterProcess(ionIoni, particle);
218
219 // muons, hadrons ions
221
222 // extra configuration
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fMinimal
#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
static G4EmParameters * Instance()
G4double MscEnergyLimit() const
void SetStepFunction(G4double v1, G4double v2)
void SetLateralDisplacement(G4bool val)
void SetVerbose(G4int val)
void SetApplyCuts(G4bool val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMscRangeFactor(G4double val)
G4EmStandardPhysics_option2(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 SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, size_t index=0)
const G4String & GetPhysicsName() const