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

#include <G4EmLivermorePhysics.hh>

+ Inheritance diagram for G4EmLivermorePhysics:

Public Member Functions

 G4EmLivermorePhysics (G4int ver=1, const G4String &name="G4EmLivermore")
 
 ~G4EmLivermorePhysics () 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 ()
 
virtual void ConstructParticle ()=0
 
virtual void ConstructProcess ()=0
 
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 35 of file G4EmLivermorePhysics.hh.

Constructor & Destructor Documentation

◆ G4EmLivermorePhysics()

G4EmLivermorePhysics::G4EmLivermorePhysics ( G4int  ver = 1,
const G4String name = "G4EmLivermore" 
)
explicit

Definition at line 111 of file G4EmLivermorePhysics.cc.

112 : G4VPhysicsConstructor(pname)
113{
114 SetVerboseLevel(ver);
116 param->SetDefaults();
117 param->SetVerbose(ver);
118 param->SetMinEnergy(100*CLHEP::eV);
119 param->SetLowestElectronEnergy(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);
132 param->SetUseICRU90Data(true);
134 param->SetMaxNIELEnergy(1*CLHEP::MeV);
136}
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseSafetyPlus
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 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 SetVerboseLevel(G4int value)

◆ ~G4EmLivermorePhysics()

G4EmLivermorePhysics::~G4EmLivermorePhysics ( )
override

Definition at line 140 of file G4EmLivermorePhysics.cc.

141{}

Member Function Documentation

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 145 of file G4EmLivermorePhysics.cc.

146{
147 // minimal set of particles for EM physics
149}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360

◆ ConstructProcess()

void G4EmLivermorePhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 153 of file G4EmLivermorePhysics.cc.

154{
155 if(verboseLevel > 1) {
156 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
157 }
161
162 // processes used by several particles
163 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
164
165 // high energy limit for e+- scattering models
166 G4double highEnergyLimit= param->MscEnergyLimit();
167 G4double livEnergyLimit = 1*CLHEP::GeV;
168
169 // nuclear stopping
170 G4double nielEnergyLimit = param->MaxNIELEnergy();
171 G4NuclearStopping* pnuc = nullptr;
172 if(nielEnergyLimit > 0.0) {
173 pnuc = new G4NuclearStopping();
174 pnuc->SetMaxKinEnergy(nielEnergyLimit);
175 }
176
177 // Add Livermore EM Processes
178
179 // Add gamma EM Processes
181 G4bool polar = param->EnablePolarisation();
182
183 // photoelectric effect - Livermore model only
186 pe->SetEmModel(peModel);
187 if(polar) {
189 }
190
191 // Compton scattering - Livermore model only
194 G4VEmModel* cModel = nullptr;
195 if(polar) {
197 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
198 } else {
199 cModel = new G4LivermoreComptonModel();
200 cModel->SetHighEnergyLimit(livEnergyLimit);
201 }
202 cs->AddEmModel(0, cModel);
203
204 // gamma conversion
207 gc->SetEmModel(convLiv);
208
209 // Rayleigh scattering
211 if(polar) {
213 }
214
215 ph->RegisterProcess(pe, particle);
216 ph->RegisterProcess(cs, particle);
217 ph->RegisterProcess(gc, particle);
218 ph->RegisterProcess(rl, particle);
219
220 // e-
221 particle = G4Electron::Electron();
222
223 // multiple and single scattering
226 msc1->SetHighEnergyLimit(highEnergyLimit);
227 msc2->SetLowEnergyLimit(highEnergyLimit);
228 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
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();
240 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
241 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
242 eioni->AddEmModel(0, theIoniLiv);
243
244 // Bremsstrahlung
250 brem->SetEmModel(br1);
251 brem->SetEmModel(br2);
252 br1->SetHighEnergyLimit(GeV);
253
255
256 // register processes
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 msc1 = new G4GoudsmitSaundersonMscModel();
267 msc2 = new G4WentzelVIModel();
268 msc1->SetHighEnergyLimit(highEnergyLimit);
269 msc2->SetLowEnergyLimit(highEnergyLimit);
270 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
271
272 ssm = new G4eCoulombScatteringModel();
273 ss = new G4CoulombScattering();
274 ss->SetEmModel(ssm);
275 ss->SetMinKinEnergy(highEnergyLimit);
276 ssm->SetLowEnergyLimit(highEnergyLimit);
277 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
278
279 // ionisation
280 eioni = new G4eIonisation();
281
282 // Bremsstrahlung from standard
283 brem = new G4eBremsstrahlung();
284 br1 = new G4SeltzerBergerModel();
285 br2 = new G4eBremsstrahlungRelModel();
288 brem->SetEmModel(br1);
289 brem->SetEmModel(br2);
290 br1->SetHighEnergyLimit(GeV);
291
292 // register processes
293 ph->RegisterProcess(eioni, particle);
294 ph->RegisterProcess(brem, particle);
295 ph->RegisterProcess(ee, particle);
296 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
297 ph->RegisterProcess(ss, particle);
298
299 // generic ion
300 particle = G4GenericIon::GenericIon();
301 G4ionIonisation* ionIoni = new G4ionIonisation();
303 ph->RegisterProcess(hmsc, particle);
304 ph->RegisterProcess(ionIoni, particle);
305 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
306
307 // muons, hadrons, ions
309
310 // extra configuration
312
313}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#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 ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
Definition: G4EmBuilder.cc:409
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:399
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
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 SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
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

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