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

#include <G4EmStandardPhysics_option2.hh>

+ Inheritance diagram for G4EmStandardPhysics_option2:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option2()

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

Definition at line 96 of file G4EmStandardPhysics_option2.cc.

98 : G4VPhysicsConstructor("G4EmStandard_opt2")
99{
100 SetVerboseLevel(ver);
102 param->SetDefaults();
103 param->SetVerbose(ver);
104 param->SetApplyCuts(false);
105 param->SetStepFunction(0.8, 1*CLHEP::mm);
106 param->SetMscRangeFactor(0.2);
107 param->SetLateralDisplacement(false);
110}
@ bElectromagnetic
static G4EmParameters * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetLateralDisplacement(G4bool val)
void SetVerbose(G4int val)
void SetApplyCuts(G4bool val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMscRangeFactor(G4double val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics_option2()

G4EmStandardPhysics_option2::~G4EmStandardPhysics_option2 ( )
override

Definition at line 114 of file G4EmStandardPhysics_option2.cc.

115{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option2::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 119 of file G4EmStandardPhysics_option2.cc.

120{
121 // minimal set of particles for EM physics
123}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysics_option2::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 127 of file G4EmStandardPhysics_option2.cc.

128{
129 if(verboseLevel > 1) {
130 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
131 }
133
135
136 // processes used by several particles
137 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
138 G4NuclearStopping* pnuc(nullptr);
139
140 // high energy limit for e+- scattering models and bremsstrahlung
141 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
142
143 // Add gamma EM Processes
145
148
149 if(G4EmParameters::Instance()->GeneralProcessActive()) {
151 sp->AddEmProcess(pee);
152 sp->AddEmProcess(new G4ComptonScattering());
153 sp->AddEmProcess(new G4GammaConversion());
155 ph->RegisterProcess(sp, particle);
156
157 } else {
158 ph->RegisterProcess(pee, particle);
159 ph->RegisterProcess(new G4ComptonScattering(), particle);
160 ph->RegisterProcess(new G4GammaConversion(), particle);
161 }
162
163 // e-
164 particle = G4Electron::Electron();
165
166 G4eIonisation* eioni = new G4eIonisation();
167
168 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
170 msc1->SetHighEnergyLimit(highEnergyLimit);
171 msc2->SetLowEnergyLimit(highEnergyLimit);
172 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
173
176 ss->SetEmModel(ssm);
177 ss->SetMinKinEnergy(highEnergyLimit);
178 ssm->SetLowEnergyLimit(highEnergyLimit);
179 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
180
181 ph->RegisterProcess(eioni, particle);
182 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
183 ph->RegisterProcess(ss, particle);
184
185 // e+
186 particle = G4Positron::Positron();
187 eioni = new G4eIonisation();
188
189 msc1 = new G4UrbanMscModel();
190 msc2 = new G4WentzelVIModel();
191 msc1->SetHighEnergyLimit(highEnergyLimit);
192 msc2->SetLowEnergyLimit(highEnergyLimit);
193 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
194
195 ssm = new G4eCoulombScatteringModel();
196 ss = new G4CoulombScattering();
197 ss->SetEmModel(ssm);
198 ss->SetMinKinEnergy(highEnergyLimit);
199 ssm->SetLowEnergyLimit(highEnergyLimit);
200 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
201
202 ph->RegisterProcess(eioni, particle);
203 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
204 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
205 ph->RegisterProcess(ss, particle);
206
207 // generic ion
208 particle = G4GenericIon::GenericIon();
209 G4ionIonisation* ionIoni = new G4ionIonisation();
210 ph->RegisterProcess(hmsc, particle);
211 ph->RegisterProcess(ionIoni, particle);
212
213 // muons, hadrons ions
215
216 // extra configuration
218}
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 ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
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)
const G4String & GetPhysicsName() const

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