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

#include <G4EmStandardPhysicsWVI.hh>

+ Inheritance diagram for G4EmStandardPhysicsWVI:

Public Member Functions

 G4EmStandardPhysicsWVI (G4int ver=1)
 
 ~G4EmStandardPhysicsWVI () 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 52 of file G4EmStandardPhysicsWVI.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysicsWVI()

G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI ( G4int ver = 1)
explicit

Definition at line 96 of file G4EmStandardPhysicsWVI.cc.

97 : G4VPhysicsConstructor("G4EmStandardWVI")
98{
100 param->SetDefaults();
101 param->SetVerbose(ver);
102 param->SetMinEnergy(10*CLHEP::eV);
103 param->SetLowestElectronEnergy(100*CLHEP::eV);
104 param->SetNumberOfBinsPerDecade(20);
106 param->SetStepFunction(0.2, 100*CLHEP::um);
107 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
108 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
109 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
110 param->SetUseMottCorrection(true);
111 param->SetMuHadLateralDisplacement(true);
112 param->SetUseICRU90Data(true);
113 param->SetMscThetaLimit(0.15);
114 param->SetFluo(true);
115 param->SetMaxNIELEnergy(1*CLHEP::MeV);
117}
@ bElectromagnetic
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetMscThetaLimit(G4double val)
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 SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
G4VPhysicsConstructor(const G4String &="")

◆ ~G4EmStandardPhysicsWVI()

G4EmStandardPhysicsWVI::~G4EmStandardPhysicsWVI ( )
override

Definition at line 121 of file G4EmStandardPhysicsWVI.cc.

122{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysicsWVI::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 126 of file G4EmStandardPhysicsWVI.cc.

127{
128 // minimal set of particles for EM physics
130}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysicsWVI::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 134 of file G4EmStandardPhysicsWVI.cc.

135{
136 if(verboseLevel > 1) {
137 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
138 }
142
143 // common processes
144 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
145
146 // nuclear stopping is enabled if th eenergy limit above zero
147 G4double nielEnergyLimit = param->MaxNIELEnergy();
148 G4NuclearStopping* pnuc = nullptr;
149 if(nielEnergyLimit > 0.0) {
150 pnuc = new G4NuclearStopping();
151 pnuc->SetMaxKinEnergy(nielEnergyLimit);
152 }
153
154 // high energy limit for e+- scattering models
155 G4double highEnergyLimit = 1*CLHEP::MeV;
156
157 // Add gamma EM processes
159
162
165
167 if(param->EnablePolarisation()) {
169 }
170
171 ph->RegisterProcess(pee, particle);
172 ph->RegisterProcess(cs, particle);
173 ph->RegisterProcess(gc, particle);
174 ph->RegisterProcess(new G4RayleighScattering(), particle);
175
176 // e-
177 particle = G4Electron::Electron();
178
180 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
182 msc1->SetHighEnergyLimit(highEnergyLimit);
183 msc2->SetLowEnergyLimit(highEnergyLimit);
184 msc->SetEmModel(msc1);
185 msc->SetEmModel(msc2);
186
189 ss->SetEmModel(ssm);
190 ss->SetMinKinEnergy(highEnergyLimit);
191 ssm->SetLowEnergyLimit(highEnergyLimit);
192 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
193
194 ph->RegisterProcess(msc, particle);
195 ph->RegisterProcess(new G4eIonisation(), particle);
196 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
197 ph->RegisterProcess(ss, particle);
198
199 // e+
200 particle = G4Positron::Positron();
201
202 msc = new G4eMultipleScattering;
203 msc1 = new G4UrbanMscModel();
204 msc2 = new G4WentzelVIModel();
205 msc1->SetHighEnergyLimit(highEnergyLimit);
206 msc2->SetLowEnergyLimit(highEnergyLimit);
207 msc->SetEmModel(msc1);
208 msc->SetEmModel(msc2);
209
210 ssm = new G4eCoulombScatteringModel();
211 ss = new G4CoulombScattering();
212 ss->SetEmModel(ssm);
213 ss->SetMinKinEnergy(highEnergyLimit);
214 ssm->SetLowEnergyLimit(highEnergyLimit);
215 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
216
219
220 ph->RegisterProcess(msc, particle);
221 ph->RegisterProcess(new G4eIonisation(), particle);
222 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
223 ph->RegisterProcess(ann, particle);
224 ph->RegisterProcess(ss, particle);
225
226 // generic ion
227 particle = G4GenericIon::GenericIon();
228 G4ionIonisation* ionIoni = new G4ionIonisation();
229 auto fluc = new G4IonFluctuations();
230 ionIoni->SetFluctModel(fluc);
232 ph->RegisterProcess(hmsc, particle);
233 ph->RegisterProcess(ionIoni, particle);
234 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
235
236 // muons, hadrons, ions
238
239 // extra configuration
241}
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()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() 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 SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
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: