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

#include <G4EmStandardPhysics.hh>

+ Inheritance diagram for G4EmStandardPhysics:

Public Member Functions

 G4EmStandardPhysics (G4int ver=1, const G4String &name="")
 
 ~G4EmStandardPhysics () 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 51 of file G4EmStandardPhysics.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics()

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

Definition at line 85 of file G4EmStandardPhysics.cc.

86 : G4VPhysicsConstructor("G4EmStandard")
87{
88 SetVerboseLevel(ver);
90 param->SetDefaults();
91 param->SetVerbose(ver);
92 param->SetGeneralProcessActive(true);
95}
@ bElectromagnetic
@ fUrbanFluctuation
static G4EmParameters * Instance()
void SetGeneralProcessActive(G4bool val)
void SetFluctuationType(G4EmFluctuationType val)
void SetVerbose(G4int val)
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics()

G4EmStandardPhysics::~G4EmStandardPhysics ( )
override

Definition at line 99 of file G4EmStandardPhysics.cc.

100{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 104 of file G4EmStandardPhysics.cc.

105{
106 // minimal set of particles for EM physics
108}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360

◆ ConstructProcess()

void G4EmStandardPhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 112 of file G4EmStandardPhysics.cc.

113{
114 if(verboseLevel > 1) {
115 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
116 }
118
121
122 // processes used by several particles
123 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
124
125 // nuclear stopping is enabled if th eenergy limit above zero
126 G4double nielEnergyLimit = param->MaxNIELEnergy();
127 G4NuclearStopping* pnuc = nullptr;
128 if(nielEnergyLimit > 0.0) {
129 pnuc = new G4NuclearStopping();
130 pnuc->SetMaxKinEnergy(nielEnergyLimit);
131 }
132
133 // high energy limit for e+- scattering models and bremsstrahlung
134 G4double highEnergyLimit = param->MscEnergyLimit();
135
136 // Add gamma EM Processes
138 G4bool polar = param->EnablePolarisation();
139
140 // Photoelectric
143 pe->SetEmModel(peModel);
144 if(polar) {
146 }
147
148 // Compton scattering
150 if(polar) {
152 }
153
154 // default Rayleigh scattering is Livermore
156 if(polar) {
158 }
159
160 if(G4EmParameters::Instance()->GeneralProcessActive()) {
162 sp->AddEmProcess(pe);
163 sp->AddEmProcess(cs);
164 sp->AddEmProcess(new G4GammaConversion());
165 sp->AddEmProcess(rl);
167 ph->RegisterProcess(sp, particle);
168
169 } else {
170 ph->RegisterProcess(pe, particle);
171 ph->RegisterProcess(cs, particle);
172 ph->RegisterProcess(new G4GammaConversion(), particle);
173 ph->RegisterProcess(rl, particle);
174 }
175
176 // e-
177 particle = G4Electron::Electron();
178
179 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
181 msc1->SetHighEnergyLimit(highEnergyLimit);
182 msc2->SetLowEnergyLimit(highEnergyLimit);
183 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
184
187 ss->SetEmModel(ssm);
188 ss->SetMinKinEnergy(highEnergyLimit);
189 ssm->SetLowEnergyLimit(highEnergyLimit);
190 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
191
192 ph->RegisterProcess(new G4eIonisation(), particle);
193 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
194 ph->RegisterProcess(ss, particle);
195
196 // e+
197 particle = G4Positron::Positron();
198
199 msc1 = new G4UrbanMscModel();
200 msc2 = new G4WentzelVIModel();
201 msc1->SetHighEnergyLimit(highEnergyLimit);
202 msc2->SetLowEnergyLimit(highEnergyLimit);
203 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
204
205 ssm = new G4eCoulombScatteringModel();
206 ss = new G4CoulombScattering();
207 ss->SetEmModel(ssm);
208 ss->SetMinKinEnergy(highEnergyLimit);
209 ssm->SetLowEnergyLimit(highEnergyLimit);
210 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
211
212 ph->RegisterProcess(new G4eIonisation(), particle);
213 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
214 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
215 ph->RegisterProcess(ss, particle);
216
217 // generic ion
218 particle = G4GenericIon::GenericIon();
219 G4ionIonisation* ionIoni = new G4ionIonisation();
220 ph->RegisterProcess(hmsc, particle);
221 ph->RegisterProcess(ionIoni, particle);
222 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
223
224 // muons, hadrons ions
226
227 // extra configuration
229}
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 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: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 SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
const G4String & GetPhysicsName() const

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