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

#include <G4EmStandardPhysics_option3.hh>

+ Inheritance diagram for G4EmStandardPhysics_option3:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option3()

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

Definition at line 100 of file G4EmStandardPhysics_option3.cc.

102 : G4VPhysicsConstructor("G4EmStandard_opt3")
103{
104 SetVerboseLevel(ver);
106 param->SetDefaults();
107 param->SetVerbose(ver);
108 param->SetGeneralProcessActive(true);
109 param->SetMinEnergy(10*CLHEP::eV);
110 param->SetLowestElectronEnergy(100*CLHEP::eV);
111 param->SetNumberOfBinsPerDecade(20);
113 param->SetUseMottCorrection(true);
114 param->SetStepFunction(0.2, 100*CLHEP::um);
115 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
116 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
117 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
119 param->SetMuHadLateralDisplacement(true);
120 param->SetLateralDisplacementAlg96(true);
121 param->SetUseICRU90Data(true);
123 param->SetFluo(true);
124 param->SetMaxNIELEnergy(1*CLHEP::MeV);
126}
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseDistanceToBoundary
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetGeneralProcessActive(G4bool val)
void SetLateralDisplacementAlg96(G4bool val)
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 SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics_option3()

G4EmStandardPhysics_option3::~G4EmStandardPhysics_option3 ( )
override

Definition at line 130 of file G4EmStandardPhysics_option3.cc.

131{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option3::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 135 of file G4EmStandardPhysics_option3.cc.

136{
137 // minimal set of particles for EM physics
139}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysics_option3::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 143 of file G4EmStandardPhysics_option3.cc.

144{
145 if(verboseLevel > 1) {
146 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
147 }
149
152
153 // processes used by several particles
154 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
155
156 // nuclear stopping is enabled if th eenergy limit above zero
157 G4double nielEnergyLimit = param->MaxNIELEnergy();
158 G4NuclearStopping* pnuc = nullptr;
159 if(nielEnergyLimit > 0.0) {
160 pnuc = new G4NuclearStopping();
161 pnuc->SetMaxKinEnergy(nielEnergyLimit);
162 }
163
164 // Add gamma EM Processes
166
169 pe->SetEmModel(peModel);
170 if(param->EnablePolarisation()) {
172 }
173
176
178 if(param->EnablePolarisation()) {
180 }
181
183 if(param->EnablePolarisation()) {
185 }
186
187 if(G4EmParameters::Instance()->GeneralProcessActive()) {
189 sp->AddEmProcess(pe);
190 sp->AddEmProcess(cs);
191 sp->AddEmProcess(gc);
192 sp->AddEmProcess(rl);
194 ph->RegisterProcess(sp, particle);
195 } else {
196 ph->RegisterProcess(pe, particle);
197 ph->RegisterProcess(cs, particle);
198 ph->RegisterProcess(gc, particle);
199 ph->RegisterProcess(rl, particle);
200 }
201
202 // e-
203 particle = G4Electron::Electron();
204
205 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
206 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
207
208 G4eIonisation* eIoni = new G4eIonisation();
209
215 brem->SetEmModel(br1);
216 brem->SetEmModel(br2);
217 br2->SetLowEnergyLimit(CLHEP::GeV);
218
220
221 ph->RegisterProcess(eIoni, particle);
222 ph->RegisterProcess(brem, particle);
223 ph->RegisterProcess(ee, particle);
224
225 // e+
226 particle = G4Positron::Positron();
227
228 msc1 = new G4UrbanMscModel();
229 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
230
231 eIoni = new G4eIonisation();
232
233 brem = new G4eBremsstrahlung();
234 br1 = new G4SeltzerBergerModel();
235 br2 = new G4eBremsstrahlungRelModel();
238 brem->SetEmModel(br1);
239 brem->SetEmModel(br2);
240 br2->SetLowEnergyLimit(CLHEP::GeV);
241
242 ph->RegisterProcess(eIoni, particle);
243 ph->RegisterProcess(brem, particle);
244 ph->RegisterProcess(ee, particle);
245 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
246
247 // generic ion
248 particle = G4GenericIon::GenericIon();
249 G4ionIonisation* ionIoni = new G4ionIonisation();
250 auto fluc = new G4IonFluctuations();
251 ionIoni->SetFluctModel(fluc);
253 ph->RegisterProcess(hmsc, particle);
254 ph->RegisterProcess(ionIoni, particle);
255 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
256
257 // muons, hadrons, ions
258 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
259
260 // extra configuration
262}
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()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() 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 SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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