Geant4 11.1.1
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 ()
 
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_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->SetMscRangeFactor(0.03);
120 param->SetMuHadLateralDisplacement(true);
121 param->SetLateralDisplacementAlg96(true);
122 param->SetUseICRU90Data(true);
124 param->SetFluo(true);
125 param->SetMaxNIELEnergy(1*CLHEP::MeV);
127}
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseSafetyPlus
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)
void SetMscRangeFactor(G4double val)
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics_option3()

G4EmStandardPhysics_option3::~G4EmStandardPhysics_option3 ( )
override

Definition at line 131 of file G4EmStandardPhysics_option3.cc.

132{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option3::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 136 of file G4EmStandardPhysics_option3.cc.

137{
138 // minimal set of particles for EM physics
140}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360

◆ ConstructProcess()

void G4EmStandardPhysics_option3::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 144 of file G4EmStandardPhysics_option3.cc.

145{
146 if(verboseLevel > 1) {
147 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
148 }
150
153
154 // processes used by several particles
155 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
156
157 // nuclear stopping is enabled if th eenergy limit above zero
158 G4double nielEnergyLimit = param->MaxNIELEnergy();
159 G4NuclearStopping* pnuc = nullptr;
160 if(nielEnergyLimit > 0.0) {
161 pnuc = new G4NuclearStopping();
162 pnuc->SetMaxKinEnergy(nielEnergyLimit);
163 }
164
165 // Add gamma EM Processes
167
170 pe->SetEmModel(peModel);
171 if(param->EnablePolarisation()) {
173 }
174
177
179 if(param->EnablePolarisation()) {
181 }
182
184 if(param->EnablePolarisation()) {
186 }
187
188 if(G4EmParameters::Instance()->GeneralProcessActive()) {
190 sp->AddEmProcess(pe);
191 sp->AddEmProcess(cs);
192 sp->AddEmProcess(gc);
193 sp->AddEmProcess(rl);
195 ph->RegisterProcess(sp, particle);
196 } else {
197 ph->RegisterProcess(pe, particle);
198 ph->RegisterProcess(cs, particle);
199 ph->RegisterProcess(gc, particle);
200 ph->RegisterProcess(rl, particle);
201 }
202
203 // e-
204 particle = G4Electron::Electron();
205
206 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
207 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
208
209 G4eIonisation* eIoni = new G4eIonisation();
210
216 brem->SetEmModel(br1);
217 brem->SetEmModel(br2);
218 br2->SetLowEnergyLimit(CLHEP::GeV);
219
221
222 ph->RegisterProcess(eIoni, particle);
223 ph->RegisterProcess(brem, particle);
224 ph->RegisterProcess(ee, particle);
225
226 // e+
227 particle = G4Positron::Positron();
228
229 msc1 = new G4UrbanMscModel();
230 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
231
232 eIoni = new G4eIonisation();
233
234 brem = new G4eBremsstrahlung();
235 br1 = new G4SeltzerBergerModel();
236 br2 = new G4eBremsstrahlungRelModel();
239 brem->SetEmModel(br1);
240 brem->SetEmModel(br2);
241 br2->SetLowEnergyLimit(CLHEP::GeV);
242
243 ph->RegisterProcess(eIoni, particle);
244 ph->RegisterProcess(brem, particle);
245 ph->RegisterProcess(ee, particle);
246 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
247
248 // generic ion
249 particle = G4GenericIon::GenericIon();
250 G4ionIonisation* ionIoni = new G4ionIonisation();
251 auto fluc = new G4IonFluctuations();
252 ionIoni->SetFluctModel(fluc);
254 ph->RegisterProcess(hmsc, particle);
255 ph->RegisterProcess(ionIoni, particle);
256 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
257
258 // muons, hadrons, ions
259 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
260
261 // extra configuration
263}
double G4double
Definition: G4Types.hh:83
#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
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 SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
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: