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

#include <G4EmStandardPhysics_option4.hh>

+ Inheritance diagram for G4EmStandardPhysics_option4:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option4()

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

Definition at line 107 of file G4EmStandardPhysics_option4.cc.

109 : G4VPhysicsConstructor("G4EmStandard_opt4")
110{
111 SetVerboseLevel(ver);
113 param->SetDefaults();
114 param->SetVerbose(ver);
115 param->SetGeneralProcessActive(true);
116 param->SetMinEnergy(100*CLHEP::eV);
117 param->SetLowestElectronEnergy(100*CLHEP::eV);
118 param->SetNumberOfBinsPerDecade(20);
120 param->SetStepFunction(0.2, 10*CLHEP::um);
121 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
122 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
123 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
124 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
125 param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
126 param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
127 param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
128 param->SetMuHadLateralDisplacement(true);
129 param->SetFluo(true);
130 param->SetUseICRU90Data(true);
132 param->SetMaxNIELEnergy(1*CLHEP::MeV);
134}
@ 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 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 SetMscSkin(G4double 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_option4()

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
override

Definition at line 138 of file G4EmStandardPhysics_option4.cc.

139{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option4::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 143 of file G4EmStandardPhysics_option4.cc.

144{
145 // minimal set of particles for EM physics
147}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360

◆ ConstructProcess()

void G4EmStandardPhysics_option4::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 151 of file G4EmStandardPhysics_option4.cc.

152{
153 if(verboseLevel > 1) {
154 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
155 }
157
160
161 // processes used by several particles
162 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
163
164 // nuclear stopping is enabled if the energy limit above zero
165 G4double nielEnergyLimit = param->MaxNIELEnergy();
166 G4NuclearStopping* pnuc = nullptr;
167 if(nielEnergyLimit > 0.0) {
168 pnuc = new G4NuclearStopping();
169 pnuc->SetMaxKinEnergy(nielEnergyLimit);
170 }
171
172 // high energy limit for e+- scattering models and bremsstrahlung
173 G4double highEnergyLimit = param->MscEnergyLimit();
174
175 // Add gamma EM Processes
177 G4bool polar = param->EnablePolarisation();
178
179 // Photoelectric
182 pe->SetEmModel(peModel);
183 if(polar) {
185 }
186
187 // Compton scattering
190 G4VEmModel* cModel = nullptr;
191 if(polar) {
192 cModel = new G4LowEPPolarizedComptonModel();
193 } else {
194 cModel = new G4LowEPComptonModel();
195 }
196 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
197 cs->AddEmModel(0, cModel);
198
199 // Gamma conversion
201 G4VEmModel* conv = new G4BetheHeitler5DModel();
202 gc->SetEmModel(conv);
203
204 // default Rayleigh scattering is Livermore
206 if(polar) {
208 }
209
210 if(param->GeneralProcessActive()) {
212 sp->AddEmProcess(pe);
213 sp->AddEmProcess(cs);
214 sp->AddEmProcess(gc);
215 sp->AddEmProcess(rl);
217 ph->RegisterProcess(sp, particle);
218 } else {
219 ph->RegisterProcess(pe, particle);
220 ph->RegisterProcess(cs, particle);
221 ph->RegisterProcess(gc, particle);
222 ph->RegisterProcess(rl, particle);
223 }
224
225 // e-
226 particle = G4Electron::Electron();
227
228 // e-/e+ msc gs with Mott-correction
229 // (Mott-correction is set through G4EmParameters)
232 msc1->SetHighEnergyLimit(highEnergyLimit);
233 msc2->SetLowEnergyLimit(highEnergyLimit);
234 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
235
238 ss->SetEmModel(ssm);
239 ss->SetMinKinEnergy(highEnergyLimit);
240 ssm->SetLowEnergyLimit(highEnergyLimit);
241 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
242
243 // ionisation
244 G4eIonisation* eioni = new G4eIonisation();
246 G4VEmModel* theIoniMod = new G4PenelopeIonisationModel();
247 theIoniMod->SetHighEnergyLimit(0.1*CLHEP::MeV);
248 eioni->AddEmModel(0, theIoniMod);
249
250 // bremsstrahlung
256 brem->SetEmModel(br1);
257 brem->SetEmModel(br2);
258 br1->SetHighEnergyLimit(CLHEP::GeV);
259
261
262 // register processes
263 ph->RegisterProcess(eioni, particle);
264 ph->RegisterProcess(brem, particle);
265 ph->RegisterProcess(ee, particle);
266 ph->RegisterProcess(ss, particle);
267
268 // e+
269 particle = G4Positron::Positron();
270
271 // e-/e+ msc gs with Mott-correction
272 // (Mott-correction is set through G4EmParameters)
273 msc1 = new G4GoudsmitSaundersonMscModel();
274 msc2 = new G4WentzelVIModel();
275 msc1->SetHighEnergyLimit(highEnergyLimit);
276 msc2->SetLowEnergyLimit(highEnergyLimit);
277 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
278
279 ssm = new G4eCoulombScatteringModel();
280 ss = new G4CoulombScattering();
281 ss->SetEmModel(ssm);
282 ss->SetMinKinEnergy(highEnergyLimit);
283 ssm->SetLowEnergyLimit(highEnergyLimit);
284 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
285
286 // ionisation
287 eioni = new G4eIonisation();
290 pen->SetHighEnergyLimit(0.1*CLHEP::MeV);
291 eioni->AddEmModel(0, pen);
292
293 // bremsstrahlung
294 brem = new G4eBremsstrahlung();
295 br1 = new G4SeltzerBergerModel();
296 br2 = new G4eBremsstrahlungRelModel();
299 brem->SetEmModel(br1);
300 brem->SetEmModel(br2);
301 br1->SetHighEnergyLimit(CLHEP::GeV);
302
303 // register processes
304 ph->RegisterProcess(eioni, particle);
305 ph->RegisterProcess(brem, particle);
306 ph->RegisterProcess(ee, particle);
307 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
308 ph->RegisterProcess(ss, particle);
309
310 // generic ion
311 particle = G4GenericIon::GenericIon();
312 G4ionIonisation* ionIoni = new G4ionIonisation();
314 ph->RegisterProcess(hmsc, particle);
315 ph->RegisterProcess(ionIoni, particle);
316 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
317
318 // muons, hadrons, ions
320
321 // extra configuration
323}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#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
G4double MscEnergyLimit() const
G4bool GeneralProcessActive() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
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 SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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