Geant4 11.3.0
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{
99 G4EmParameters* param = G4EmParameters::Instance();
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 ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysicsWVI::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 125 of file G4EmStandardPhysicsWVI.cc.

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

◆ ConstructProcess()

void G4EmStandardPhysicsWVI::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 133 of file G4EmStandardPhysicsWVI.cc.

134{
135 if(verboseLevel > 1) {
136 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
137 }
139 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
140 G4EmParameters* param = G4EmParameters::Instance();
141
142 // common processes
143 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
144
145 // nuclear stopping is enabled if th eenergy limit above zero
146 G4double nielEnergyLimit = param->MaxNIELEnergy();
147 G4NuclearStopping* pnuc = nullptr;
148 if(nielEnergyLimit > 0.0) {
149 pnuc = new G4NuclearStopping();
150 pnuc->SetMaxKinEnergy(nielEnergyLimit);
151 }
152
153 // high energy limit for e+- scattering models
154 G4double highEnergyLimit = 1*CLHEP::MeV;
155
156 // Add gamma EM processes
157 G4ParticleDefinition* particle = G4Gamma::Gamma();
158
159 G4PhotoElectricEffect* pee = new G4PhotoElectricEffect();
160 pee->SetEmModel(new G4LivermorePhotoElectricModel());
161
162 G4ComptonScattering* cs = new G4ComptonScattering;
163 cs->SetEmModel(new G4KleinNishinaModel());
164
165 G4GammaConversion* gc = new G4GammaConversion();
166 if(param->EnablePolarisation()) {
167 gc->SetEmModel(new G4BetheHeitler5DModel());
168 }
169
170 ph->RegisterProcess(pee, particle);
171 ph->RegisterProcess(cs, particle);
172 ph->RegisterProcess(gc, particle);
173 ph->RegisterProcess(new G4RayleighScattering(), particle);
174
175 // e-
176 particle = G4Electron::Electron();
177
178 G4eMultipleScattering* msc = new G4eMultipleScattering;
179 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
180 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
181 msc1->SetHighEnergyLimit(highEnergyLimit);
182 msc2->SetLowEnergyLimit(highEnergyLimit);
183 msc->SetEmModel(msc1);
184 msc->SetEmModel(msc2);
185
186 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
187 G4CoulombScattering* ss = new G4CoulombScattering();
188 ss->SetEmModel(ssm);
189 ss->SetMinKinEnergy(highEnergyLimit);
190 ssm->SetLowEnergyLimit(highEnergyLimit);
191 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
192
193 ph->RegisterProcess(msc, particle);
194 ph->RegisterProcess(new G4eIonisation(), particle);
195 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
196 ph->RegisterProcess(ss, particle);
197
198 // e+
199 particle = G4Positron::Positron();
200
201 msc = new G4eMultipleScattering;
202 msc1 = new G4UrbanMscModel();
203 msc2 = new G4WentzelVIModel();
204 msc1->SetHighEnergyLimit(highEnergyLimit);
205 msc2->SetLowEnergyLimit(highEnergyLimit);
206 msc->SetEmModel(msc1);
207 msc->SetEmModel(msc2);
208
209 ssm = new G4eCoulombScatteringModel();
210 ss = new G4CoulombScattering();
211 ss->SetEmModel(ssm);
212 ss->SetMinKinEnergy(highEnergyLimit);
213 ssm->SetLowEnergyLimit(highEnergyLimit);
214 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
215
216 G4eplusAnnihilation* ann = new G4eplusAnnihilation();
217 ann->SetEmModel(new G4eplusTo2or3GammaModel());
218
219 ph->RegisterProcess(msc, particle);
220 ph->RegisterProcess(new G4eIonisation(), particle);
221 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
222 ph->RegisterProcess(ann, particle);
223 ph->RegisterProcess(ss, particle);
224
225 // generic ion
226 particle = G4GenericIon::GenericIon();
227 G4ionIonisation* ionIoni = new G4ionIonisation();
228 auto fluc = new G4AtimaFluctuations();
229 ionIoni->SetFluctModel(fluc);
230 ionIoni->SetEmModel(new G4AtimaEnergyLossModel());
231 ph->RegisterProcess(hmsc, particle);
232 ph->RegisterProcess(ionIoni, particle);
233 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
234
235 // muons, hadrons, ions
237
238 // extra configuration
239 G4EmModelActivator mact(GetPhysicsName());
240}
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: