Geant4 10.7.0
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
 
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 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 95 of file G4EmStandardPhysics_option3.cc.

97 : G4VPhysicsConstructor("G4EmStandard_opt3"), verbose(ver)
98{
100 param->SetDefaults();
101 param->SetVerbose(verbose);
102 param->SetMinEnergy(10*CLHEP::eV);
103 param->SetLowestElectronEnergy(100*CLHEP::eV);
104 param->SetNumberOfBinsPerDecade(20);
106 param->SetUseMottCorrection(true);
107 param->SetStepFunction(0.2, 100*CLHEP::um);
108 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
109 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
110 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
112 param->SetMuHadLateralDisplacement(true);
113 param->SetLateralDisplacementAlg96(true);
114 param->SetUseICRU90Data(true);
115 param->SetFluo(true);
116 param->SetMaxNIELEnergy(1*CLHEP::MeV);
118}
@ bElectromagnetic
@ fUseDistanceToBoundary
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetLateralDisplacementAlg96(G4bool 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 SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)

◆ ~G4EmStandardPhysics_option3()

G4EmStandardPhysics_option3::~G4EmStandardPhysics_option3 ( )
override

Definition at line 122 of file G4EmStandardPhysics_option3.cc.

123{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option3::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 127 of file G4EmStandardPhysics_option3.cc.

128{
129 // minimal set of particles for EM physics
131}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234

◆ ConstructProcess()

void G4EmStandardPhysics_option3::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 135 of file G4EmStandardPhysics_option3.cc.

136{
137 if(verbose > 1) {
138 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
139 }
141
143
144 // processes used by several particles
146 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
147
148 // nuclear stopping
149 G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
151 pnuc->SetMaxKinEnergy(nielEnergyLimit);
152
153 // Add gamma EM Processes
155
158
161
162 if(G4EmParameters::Instance()->GeneralProcessActive()) {
164 sp->AddEmProcess(pee);
165 sp->AddEmProcess(cs);
166 sp->AddEmProcess(new G4GammaConversion());
167 sp->AddEmProcess(new G4RayleighScattering());
169 ph->RegisterProcess(sp, particle);
170 } else {
171 ph->RegisterProcess(pee,particle);
172 ph->RegisterProcess(cs, particle);
173 ph->RegisterProcess(new G4GammaConversion(), particle);
174 ph->RegisterProcess(new G4RayleighScattering(), particle);
175 }
176
177 // e-
178 particle = G4Electron::Electron();
179
181 G4eIonisation* eIoni = new G4eIonisation();
182
188 brem->SetEmModel(br1);
189 brem->SetEmModel(br2);
190 br2->SetLowEnergyLimit(CLHEP::GeV);
191
192 ph->RegisterProcess(msc, particle);
193 ph->RegisterProcess(eIoni, particle);
194 ph->RegisterProcess(brem, particle);
195 ph->RegisterProcess(ee, particle);
196
197 // e+
198 particle = G4Positron::Positron();
199
200 msc = new G4eMultipleScattering();
201 eIoni = new G4eIonisation();
202
203 brem = new G4eBremsstrahlung();
204 br1 = new G4SeltzerBergerModel();
205 br2 = new G4eBremsstrahlungRelModel();
208 brem->SetEmModel(br1);
209 brem->SetEmModel(br2);
210 br2->SetLowEnergyLimit(CLHEP::GeV);
211
212 ph->RegisterProcess(msc, particle);
213 ph->RegisterProcess(eIoni, particle);
214 ph->RegisterProcess(brem, particle);
215 ph->RegisterProcess(ee, particle);
216 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
217
218 // generic ion
219 particle = G4GenericIon::GenericIon();
220 G4ionIonisation* ionIoni = new G4ionIonisation();
222 ph->RegisterProcess(hmsc, particle);
223 ph->RegisterProcess(ionIoni, particle);
224 ph->RegisterProcess(pnuc, particle);
225
226 // muons, hadrons, ions
227 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
228
229 // extra configuration
231}
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:170
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
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:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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