Geant4 10.7.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 ()
 
virtual void ConstructParticle ()=0
 
virtual void ConstructProcess ()=0
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 

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
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4int g4vpcInstanceID
 
- 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 95 of file G4EmStandardPhysicsWVI.cc.

96 : G4VPhysicsConstructor("G4EmStandardWVI"), verbose(ver)
97{
99 param->SetDefaults();
100 param->SetVerbose(verbose);
101 param->SetMinEnergy(10*CLHEP::eV);
102 param->SetLowestElectronEnergy(10*CLHEP::eV);
103 param->SetNumberOfBinsPerDecade(20);
105 param->SetStepFunction(0.2, 100*CLHEP::um);
106 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
107 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
108 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
109 param->SetUseMottCorrection(true);
110 param->SetMuHadLateralDisplacement(true);
111 param->SetMscThetaLimit(0.15);
112 param->SetFluo(true);
114}
@ 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 SetStepFunctionIons(G4double v1, G4double v2)
void SetUseMottCorrection(G4bool val)

◆ ~G4EmStandardPhysicsWVI()

G4EmStandardPhysicsWVI::~G4EmStandardPhysicsWVI ( )
override

Definition at line 118 of file G4EmStandardPhysicsWVI.cc.

119{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysicsWVI::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 123 of file G4EmStandardPhysicsWVI.cc.

124{
125 // minimal set of particles for EM physics
127}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234

◆ ConstructProcess()

void G4EmStandardPhysicsWVI::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 131 of file G4EmStandardPhysicsWVI.cc.

132{
133 if(verbose > 1) {
134 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
135 }
138
139 // common processes
140 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
141 G4NuclearStopping* pnuc(nullptr);
142
143 // Add gamma EM processes
145
148
151
152 ph->RegisterProcess(pee, particle);
153 ph->RegisterProcess(cs, particle);
154 ph->RegisterProcess(new G4GammaConversion(), particle);
155 ph->RegisterProcess(new G4RayleighScattering(), particle);
156
157 // e-
158 particle = G4Electron::Electron();
159
161 msc->SetEmModel(new G4WentzelVIModel());
163
164 ph->RegisterProcess(msc, particle);
165 ph->RegisterProcess(new G4eIonisation(), particle);
166 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
167 ph->RegisterProcess(ss, particle);
168
169 // e+
170 particle = G4Positron::Positron();
171
172 msc = new G4eMultipleScattering;
173 msc->SetEmModel(new G4WentzelVIModel());
174 ss = new G4CoulombScattering();
175
178
179 ph->RegisterProcess(msc, particle);
180 ph->RegisterProcess(new G4eIonisation(), particle);
181 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
182 ph->RegisterProcess(ann, particle);
183 ph->RegisterProcess(ss, particle);
184
185 // generic ion
186 particle = G4GenericIon::GenericIon();
187 G4ionIonisation* ionIoni = new G4ionIonisation();
188 ionIoni->SetEmModel(new G4BraggIonModel(),0);
189 ionIoni->SetEmModel(new G4AtimaEnergyLossModel(),1);
190 ionIoni->SetFluctModel(new G4AtimaFluctuations());
191 ph->RegisterProcess(hmsc, particle);
192 ph->RegisterProcess(ionIoni, particle);
193
194 // muons, hadrons, ions
196
197 // extra configuration
199}
#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:170
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetEmModel(G4VEmModel *, G4int index=0)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, size_t index=0)
const G4String & GetPhysicsName() const

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