Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmLowEPPhysics Class Reference

#include <G4EmLowEPPhysics.hh>

+ Inheritance diagram for G4EmLowEPPhysics:

Public Member Functions

 G4EmLowEPPhysics (G4int ver=1, const G4String &name="")
 
 ~G4EmLowEPPhysics () override
 
void ConstructParticle () override
 
void ConstructProcess () override
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
PhysicsBuilder_V GetBuilders () const
 
void AddBuilder (G4PhysicsBuilderInterface *bld)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel = 0
 
G4String namePhysics = ""
 
G4int typePhysics = 0
 
G4ParticleTabletheParticleTable = nullptr
 
G4int g4vpcInstanceID = 0
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 48 of file G4EmLowEPPhysics.hh.

Constructor & Destructor Documentation

◆ G4EmLowEPPhysics()

G4EmLowEPPhysics::G4EmLowEPPhysics ( G4int ver = 1,
const G4String & name = "" )
explicit

Definition at line 121 of file G4EmLowEPPhysics.cc.

122 : G4VPhysicsConstructor("G4EmLowEPPhysics")
123{
124 SetVerboseLevel(ver);
126 param->SetDefaults();
127 param->SetVerbose(ver);
128 param->SetMinEnergy(100*CLHEP::eV);
129 param->SetLowestElectronEnergy(100*CLHEP::eV);
130 param->SetNumberOfBinsPerDecade(20);
132 param->SetStepFunction(0.2, 100*CLHEP::um);
133 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
134 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
135 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
136 param->SetUseMottCorrection(true);
137 param->SetMscRangeFactor(0.04);
138 param->SetMuHadLateralDisplacement(true);
139 param->SetFluo(true);
140 param->SetAuger(true);
141 param->SetUseICRU90Data(true);
143}
@ bElectromagnetic
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
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 SetStepFunctionIons(G4double v1, G4double v2)
void SetAuger(G4bool val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmLowEPPhysics()

G4EmLowEPPhysics::~G4EmLowEPPhysics ( )
override

Definition at line 147 of file G4EmLowEPPhysics.cc.

148{}

Member Function Documentation

◆ ConstructParticle()

void G4EmLowEPPhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 152 of file G4EmLowEPPhysics.cc.

153{
154 // minimal set of particles for EM physics
156}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmLowEPPhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 160 of file G4EmLowEPPhysics.cc.

161{
162 if(verboseLevel > 1) {
163 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
164 }
166
169
170 // common processes
171 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
172
173 // nuclear stopping
174 G4double nielEnergyLimit = param->MaxNIELEnergy();
175 G4NuclearStopping* pnuc = nullptr;
176 if(nielEnergyLimit > 0.0) {
177 pnuc = new G4NuclearStopping();
178 pnuc->SetMaxKinEnergy(nielEnergyLimit);
179 }
180
181 // high energy limit for e+- scattering models and bremsstrahlung
182 G4double highEnergyLimit = param->MscEnergyLimit();
183
184 // Add gamma EM Processes
186
187 // Photoelectric absorption
189 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
191 pe->SetEmModel(theLivermorePEModel);
192
193 // Compton scattering - Polarised Monash model
195 G4VEmModel* theLowEPCSModel = new G4LowEPPolarizedComptonModel();
196 cs->SetEmModel(theLowEPCSModel);
197
198 // gamma conversion - 5D model below 80 GeV with Livermore x-sections
199 G4GammaConversion* theGammaConversion = new G4GammaConversion();
200 G4VEmModel* conv = new G4BetheHeitler5DModel();
201 theGammaConversion->SetEmModel(conv);
202
203 // default Rayleigh scattering is Livermore
204 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
205 G4VEmModel* theLivermorePRSModel = new G4LivermorePolarizedRayleighModel();
206 theRayleigh->SetEmModel(theLivermorePRSModel);
207
208 ph->RegisterProcess(pe, particle);
209 ph->RegisterProcess(cs, particle);
210 ph->RegisterProcess(theGammaConversion, particle);
211 ph->RegisterProcess(theRayleigh, particle);
212
213 // e-
214 particle = G4Electron::Electron();
215
216 // multiple scattering
220 msc1->SetHighEnergyLimit(highEnergyLimit);
221 msc2->SetLowEnergyLimit(highEnergyLimit);
222 msc->SetEmModel(msc1);
223 msc->SetEmModel(msc2);
224
227 ss->SetEmModel(ssm);
228 ss->SetMinKinEnergy(highEnergyLimit);
229 ssm->SetLowEnergyLimit(highEnergyLimit);
230 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
231
232 // Ionisation - Livermore should be used only for low energies
233 G4eIonisation* eioni = new G4eIonisation();
235 theIoniLivermore->SetHighEnergyLimit(0.1*CLHEP::MeV);
236 eioni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
237
238 // Bremsstrahlung
244 brem->SetEmModel(br1);
245 brem->SetEmModel(br2);
246 br1->SetHighEnergyLimit(GeV);
247
249
250 // register processes
251 ph->RegisterProcess(msc, particle);
252 ph->RegisterProcess(ss, particle);
253 ph->RegisterProcess(eioni, particle);
254 ph->RegisterProcess(brem, particle);
255 ph->RegisterProcess(ee, particle);
256
257 // e+
258 particle = G4Positron::Positron();
259
260 // multiple scattering
261 msc = new G4eMultipleScattering();
262 msc1 = new G4GoudsmitSaundersonMscModel();
263 msc2 = new G4WentzelVIModel();
264 msc1->SetHighEnergyLimit(highEnergyLimit);
265 msc2->SetLowEnergyLimit(highEnergyLimit);
266 msc->SetEmModel(msc1);
267 msc->SetEmModel(msc2);
268
269 ssm = new G4eCoulombScatteringModel();
270 ss = new G4CoulombScattering();
271 ss->SetEmModel(ssm);
272 ss->SetMinKinEnergy(highEnergyLimit);
273 ssm->SetLowEnergyLimit(highEnergyLimit);
274 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
275
276 // Standard ionisation
277 eioni = new G4eIonisation();
279 pen->SetHighEnergyLimit(0.1*MeV);
280 eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
281
282 // Bremsstrahlung
283 brem = new G4eBremsstrahlung();
284 br1 = new G4SeltzerBergerModel();
285 br2 = new G4eBremsstrahlungRelModel();
288 brem->SetEmModel(br1);
289 brem->SetEmModel(br2);
290 br1->SetHighEnergyLimit(CLHEP::GeV);
291
292 // register processes
293 ph->RegisterProcess(msc, particle);
294 ph->RegisterProcess(ss, particle);
295 ph->RegisterProcess(eioni, particle);
296 ph->RegisterProcess(brem, particle);
297 ph->RegisterProcess(ee, particle);
298 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
299
300 // generic ion
301 particle = G4GenericIon::GenericIon();
302
303 G4ionIonisation* ionIoni = new G4ionIonisation();
305 ionIoni->SetEmModel(mod);
306
307 ph->RegisterProcess(hmsc, particle);
308 ph->RegisterProcess(ionIoni, particle);
309 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
310
311 // muons, hadrons ions
313
314 // extra configuration
316}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
static void PrepareEMPhysics()
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetMinKinEnergy(G4double e)
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

The documentation for this class was generated from the following files: