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

#include <G4EmParameters.hh>

Public Member Functions

 ~G4EmParameters ()
void SetDefaults ()
void StreamInfo (std::ostream &os) const
void Dump ()
void SetLossFluctuations (G4bool val)
G4bool LossFluctuation () const
void SetBuildCSDARange (G4bool val)
G4bool BuildCSDARange () const
void SetLPM (G4bool val)
G4bool LPM () const
void SetUseCutAsFinalRange (G4bool val)
G4bool UseCutAsFinalRange () const
void SetApplyCuts (G4bool val)
G4bool ApplyCuts () const
void SetFluo (G4bool val)
G4bool Fluo () const
G4EmFluoDirectory FluoDirectory () const
void SetFluoDirectory (G4EmFluoDirectory)
void SetBeardenFluoDir (G4bool val)
void SetANSTOFluoDir (G4bool val)
void SetXDB_EADLFluoDir (G4bool val)
G4bool BeardenFluoDir ()
G4bool ANSTOFluoDir ()
void SetAuger (G4bool val)
void SetAugerCascade (G4bool val)
G4bool Auger () const
G4bool AugerCascade () const
void SetPixe (G4bool val)
G4bool Pixe () const
void SetDeexcitationIgnoreCut (G4bool val)
G4bool DeexcitationIgnoreCut () const
void SetLateralDisplacement (G4bool val)
G4bool LateralDisplacement () const
void SetLateralDisplacementAlg96 (G4bool val)
G4bool LateralDisplacementAlg96 () const
void SetMuHadLateralDisplacement (G4bool val)
G4bool MuHadLateralDisplacement () const
void ActivateAngularGeneratorForIonisation (G4bool val)
G4bool UseAngularGeneratorForIonisation () const
void SetUseMottCorrection (G4bool val)
G4bool UseMottCorrection () const
void SetIntegral (G4bool val)
G4bool Integral () const
void SetBirksActive (G4bool val)
G4bool BirksActive () const
void SetUseICRU90Data (G4bool val)
G4bool UseICRU90Data () const
void SetFluctuationType (G4EmFluctuationType val)
G4EmFluctuationType FluctuationType () const
void SetPositronAtRestModelType (G4PositronAtRestModelType val)
G4PositronAtRestModelType PositronAtRestModelType () const
void SetDNAFast (G4bool val)
G4bool DNAFast () const
void SetDNAStationary (G4bool val)
G4bool DNAStationary () const
void SetDNAElectronMsc (G4bool val)
G4bool DNAElectronMsc () const
void SetGeneralProcessActive (G4bool val)
G4bool GeneralProcessActive () const
void SetEnableSamplingTable (G4bool val)
G4bool EnableSamplingTable () const
void SetEnablePolarisation (G4bool val)
G4bool EnablePolarisation () const
G4bool GetDirectionalSplitting () const
void SetDirectionalSplitting (G4bool v)
G4bool QuantumEntanglement () const
void SetQuantumEntanglement (G4bool v)
G4bool RetrieveMuDataFromFile () const
void SetRetrieveMuDataFromFile (G4bool v)
G4bool PhotoeffectBelowKShell () const
void SetPhotoeffectBelowKShell (G4bool v)
G4bool MscPositronCorrection () const
void SetMscPositronCorrection (G4bool v)
G4bool UseEPICS2017XS () const
void SetUseEPICS2017XS (G4bool v)
G4bool Use3GammaAnnihilationOnFly () const
void Set3GammaAnnihilationOnFly (G4bool v)
G4bool UseRiGePairProductionModel () const
void SetUseRiGePairProductionModel (G4bool v)
void SetOnIsolated (G4bool val)
G4bool OnIsolated () const
void ActivateDNA ()
void SetIsPrintedFlag (G4bool val)
G4bool IsPrintLocked () const
void SetMinEnergy (G4double val)
G4double MinKinEnergy () const
void SetMaxEnergy (G4double val)
G4double MaxKinEnergy () const
void SetMaxEnergyForCSDARange (G4double val)
G4double MaxEnergyForCSDARange () const
void SetLowestElectronEnergy (G4double val)
G4double LowestElectronEnergy () const
void SetLowestMuHadEnergy (G4double val)
G4double LowestMuHadEnergy () const
void SetLowestTripletEnergy (G4double val)
G4double LowestTripletEnergy () const
void SetLinearLossLimit (G4double val)
G4double LinearLossLimit () const
void SetBremsstrahlungTh (G4double val)
G4double BremsstrahlungTh () const
void SetMuHadBremsstrahlungTh (G4double val)
G4double MuHadBremsstrahlungTh () const
void SetLambdaFactor (G4double val)
G4double LambdaFactor () const
void SetFactorForAngleLimit (G4double val)
G4double FactorForAngleLimit () const
void SetMscThetaLimit (G4double val)
G4double MscThetaLimit () const
void SetMscEnergyLimit (G4double val)
G4double MscEnergyLimit () const
void SetMscRangeFactor (G4double val)
G4double MscRangeFactor () const
void SetMscMuHadRangeFactor (G4double val)
G4double MscMuHadRangeFactor () const
void SetMscGeomFactor (G4double val)
G4double MscGeomFactor () const
void SetMscSafetyFactor (G4double val)
G4double MscSafetyFactor () const
void SetMscLambdaLimit (G4double val)
G4double MscLambdaLimit () const
void SetMscSkin (G4double val)
G4double MscSkin () const
void SetScreeningFactor (G4double val)
G4double ScreeningFactor () const
void SetMaxNIELEnergy (G4double val)
G4double MaxNIELEnergy () const
void SetMaxEnergyFor5DMuPair (G4double val)
G4double MaxEnergyFor5DMuPair () const
void SetMaxDNAElectronEnergy (G4double val)
G4double MaxDNAElectronEnergy () const
void SetMaxDNAIonEnergy (G4double val)
G4double MaxDNAIonEnergy () const
void SetStepFunction (G4double v1, G4double v2)
void SetStepFunctionMuHad (G4double v1, G4double v2)
void SetStepFunctionLightIons (G4double v1, G4double v2)
void SetStepFunctionIons (G4double v1, G4double v2)
void FillStepFunction (const G4ParticleDefinition *, G4VEnergyLossProcess *) const
void SetDirectionalSplittingRadius (G4double r)
G4double GetDirectionalSplittingRadius ()
void SetDirectionalSplittingTarget (const G4ThreeVector &v)
G4ThreeVector GetDirectionalSplittingTarget () const
void SetNumberOfBinsPerDecade (G4int val)
G4int NumberOfBinsPerDecade () const
G4int NumberOfBins () const
void SetVerbose (G4int val)
G4int Verbose () const
void SetWorkerVerbose (G4int val)
G4int WorkerVerbose () const
void SetNumberForFreeVector (G4int val)
G4int NumberForFreeVector () const
void SetTransportationWithMsc (G4TransportationWithMscType val)
G4TransportationWithMscType TransportationWithMsc () const
void SetMscStepLimitType (G4MscStepLimitType val)
G4MscStepLimitType MscStepLimitType () const
void SetMscMuHadStepLimitType (G4MscStepLimitType val)
G4MscStepLimitType MscMuHadStepLimitType () const
void SetSingleScatteringType (G4eSingleScatteringType val)
G4eSingleScatteringType SingleScatteringType () const
void SetNuclearFormfactorType (G4NuclearFormfactorType val)
G4NuclearFormfactorType NuclearFormfactorType () const
void SetDNAeSolvationSubType (G4DNAModelSubType val)
G4DNAModelSubType DNAeSolvationSubType () const
void SetTimeStepModel (const G4ChemTimeStepModel &model)
G4ChemTimeStepModel GetTimeStepModel () const
void SetConversionType (G4int val)
G4int GetConversionType () const
void SetPIXECrossSectionModel (const G4String &)
const G4StringPIXECrossSectionModel ()
void SetPIXEElectronCrossSectionModel (const G4String &)
const G4StringPIXEElectronCrossSectionModel ()
void SetLivermoreDataDir (const G4String &)
const G4StringLivermoreDataDir ()
void AddPAIModel (const G4String &particle, const G4String &region, const G4String &type)
const std::vector< G4String > & ParticlesPAI () const
const std::vector< G4String > & RegionsPAI () const
const std::vector< G4String > & TypesPAI () const
void AddMicroElec (const G4String &region)
const std::vector< G4String > & RegionsMicroElec () const
void AddDNA (const G4String &region, const G4String &type)
const std::vector< G4String > & RegionsDNA () const
const std::vector< G4String > & TypesDNA () const
void AddPhysics (const G4String &region, const G4String &type)
const std::vector< G4String > & RegionsPhysics () const
const std::vector< G4String > & TypesPhysics () const
void SetSubCutRegion (const G4String &region="")
void SetDeexActiveRegion (const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
void SetProcessBiasingFactor (const G4String &procname, G4double val, G4bool wflag)
void ActivateForcedInteraction (const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void ActivateSecondaryBiasing (const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetEmSaturation (G4EmSaturation *)
G4EmSaturationGetEmSaturation ()
void SetFluctuationsForRegion (const G4String &regionName, G4bool flag)
void DefineRegParamForLoss (G4VEnergyLossProcess *) const
void DefineRegParamForEM (G4VEmProcess *) const
void DefineRegParamForDeex (G4VAtomDeexcitation *) const
void DefineFluctuationFlags (std::vector< G4bool > *theFluctFlags)
const G4StringGetDirLEDATA () const
 G4EmParameters (G4EmParameters &)=delete
G4EmParametersoperator= (const G4EmParameters &right)=delete

Static Public Member Functions

static G4EmParametersInstance ()

Friends

std::ostream & operator<< (std::ostream &os, const G4EmParameters &par)

Detailed Description

Definition at line 103 of file G4EmParameters.hh.

Constructor & Destructor Documentation

◆ ~G4EmParameters()

G4EmParameters::~G4EmParameters ( )

Definition at line 86 of file G4EmParameters.cc.

87{
88 delete theMessenger;
89 delete fBParameters;
90 delete fCParameters;
91 delete emSaturation;
92}

◆ G4EmParameters()

G4EmParameters::G4EmParameters ( G4EmParameters & )
delete

Member Function Documentation

◆ ActivateAngularGeneratorForIonisation()

◆ ActivateDNA()

void G4EmParameters::ActivateDNA ( )

Definition at line 576 of file G4EmParameters.cc.

577{
578 if(IsLocked()) { return; }
579 fDNA = true;
580}

Referenced by AddDNA(), G4EmDNAPhysics::G4EmDNAPhysics(), SetDNAElectronMsc(), SetDNAeSolvationSubType(), SetDNAFast(), and SetDNAStationary().

◆ ActivateForcedInteraction()

void G4EmParameters::ActivateForcedInteraction ( const G4String & procname,
const G4String & region,
G4double length,
G4bool wflag )

Definition at line 1305 of file G4EmParameters.cc.

1309{
1310 if(IsLocked() && !gener) { return; }
1311 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1312}

◆ ActivateSecondaryBiasing()

void G4EmParameters::ActivateSecondaryBiasing ( const G4String & name,
const G4String & region,
G4double factor,
G4double energyLimit )

Definition at line 1315 of file G4EmParameters.cc.

1319{
1320 if(IsLocked()) { return; }
1321 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1322}

◆ AddDNA()

void G4EmParameters::AddDNA ( const G4String & region,
const G4String & type )

Definition at line 1249 of file G4EmParameters.cc.

1250{
1251 if(IsLocked()) { return; }
1252 fCParameters->AddDNA(region, type);
1253 ActivateDNA();
1254}

◆ AddMicroElec()

void G4EmParameters::AddMicroElec ( const G4String & region)

Definition at line 1238 of file G4EmParameters.cc.

1239{
1240 if(IsLocked()) { return; }
1241 fCParameters->AddMicroElec(region);
1242}

◆ AddPAIModel()

void G4EmParameters::AddPAIModel ( const G4String & particle,
const G4String & region,
const G4String & type )

Definition at line 1215 of file G4EmParameters.cc.

1218{
1219 if(IsLocked()) { return; }
1220 fBParameters->AddPAIModel(particle, region, type);
1221}

◆ AddPhysics()

void G4EmParameters::AddPhysics ( const G4String & region,
const G4String & type )

Definition at line 1266 of file G4EmParameters.cc.

1267{
1268 if(IsLocked()) { return; }
1269 fBParameters->AddPhysics(region, type);
1270}

Referenced by LBE::ConstructGeneral().

◆ ANSTOFluoDir()

G4bool G4EmParameters::ANSTOFluoDir ( )

Definition at line 299 of file G4EmParameters.cc.

300{
301 auto dir = fCParameters->FluoDirectory();
302 return (dir == fluoANSTO);
303}
@ fluoANSTO

Referenced by G4VRadioactiveDecay::StreamInfo().

◆ ApplyCuts()

G4bool G4EmParameters::ApplyCuts ( ) const

Definition at line 242 of file G4EmParameters.cc.

243{
244 return applyCuts;
245}

◆ Auger()

G4bool G4EmParameters::Auger ( ) const

Definition at line 305 of file G4EmParameters.cc.

306{
307 return fCParameters->Auger();
308}

Referenced by AugerCascade(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), and G4VRadioactiveDecay::StreamInfo().

◆ AugerCascade()

G4bool G4EmParameters::AugerCascade ( ) const
inline

Definition at line 150 of file G4EmParameters.hh.

150{ return Auger(); }
G4bool Auger() const

◆ BeardenFluoDir()

G4bool G4EmParameters::BeardenFluoDir ( )

Definition at line 293 of file G4EmParameters.cc.

294{
295 auto dir = fCParameters->FluoDirectory();
296 return (dir == fluoBearden);
297}
@ fluoBearden

Referenced by G4VRadioactiveDecay::StreamInfo().

◆ BirksActive()

G4bool G4EmParameters::BirksActive ( ) const

Definition at line 416 of file G4EmParameters.cc.

417{
418 return birks;
419}

◆ BremsstrahlungTh()

G4double G4EmParameters::BremsstrahlungTh ( ) const

Definition at line 747 of file G4EmParameters.cc.

748{
749 return bremsTh;
750}

Referenced by G4eBremsstrahlung::InitialiseEnergyLossProcess(), and G4eBremsstrahlung::StreamProcessInfo().

◆ BuildCSDARange()

G4bool G4EmParameters::BuildCSDARange ( ) const

Definition at line 209 of file G4EmParameters.cc.

210{
211 return buildCSDARange;
212}

◆ DeexcitationIgnoreCut()

G4bool G4EmParameters::DeexcitationIgnoreCut ( ) const

Definition at line 327 of file G4EmParameters.cc.

328{
329 return fCParameters->DeexcitationIgnoreCut();
330}

Referenced by G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), and G4VRadioactiveDecay::StreamInfo().

◆ DefineFluctuationFlags()

void G4EmParameters::DefineFluctuationFlags ( std::vector< G4bool > * theFluctFlags)

Definition at line 1397 of file G4EmParameters.cc.

1398{
1399 if (fluctRegions.empty()) { return; }
1400 G4EmUtility::FillFluctFlags(fluctRegions, theFlags);
1401}
static void FillFluctFlags(std::vector< std::pair< G4String, G4bool > > &reg, std::vector< G4bool > *flags)

◆ DefineRegParamForDeex()

void G4EmParameters::DefineRegParamForDeex ( G4VAtomDeexcitation * ptr) const

Definition at line 1377 of file G4EmParameters.cc.

1378{
1379 fCParameters->DefineRegParamForDeex(ptr);
1380}

Referenced by G4VAtomDeexcitation::InitialiseAtomicDeexcitation().

◆ DefineRegParamForEM()

void G4EmParameters::DefineRegParamForEM ( G4VEmProcess * ptr) const

Definition at line 1329 of file G4EmParameters.cc.

1330{
1331 fBParameters->DefineRegParamForEM(ptr);
1332}

◆ DefineRegParamForLoss()

void G4EmParameters::DefineRegParamForLoss ( G4VEnergyLossProcess * ptr) const

Definition at line 1324 of file G4EmParameters.cc.

1325{
1326 fBParameters->DefineRegParamForLoss(ptr);
1327}

◆ DNAElectronMsc()

G4bool G4EmParameters::DNAElectronMsc ( ) const

Definition at line 463 of file G4EmParameters.cc.

464{
465 return fCParameters->DNAElectronMsc();
466}

◆ DNAeSolvationSubType()

G4DNAModelSubType G4EmParameters::DNAeSolvationSubType ( ) const

Definition at line 1161 of file G4EmParameters.cc.

1162{
1163 return fCParameters->DNAeSolvationSubType();
1164}

Referenced by G4CT_COUNT_IMPL(), and G4DNASolvationModelFactory::GetMacroDefinedModel().

◆ DNAFast()

◆ DNAStationary()

◆ Dump()

void G4EmParameters::Dump ( )

Definition at line 1590 of file G4EmParameters.cc.

1591{
1592 if(fIsPrinted) return;
1593
1594#ifdef G4MULTITHREADED
1595 G4MUTEXLOCK(&emParametersMutex);
1596#endif
1598#ifdef G4MULTITHREADED
1599 G4MUTEXUNLOCK(&emParametersMutex);
1600#endif
1601}
#define G4MUTEXLOCK(mutex)
#define G4MUTEXUNLOCK(mutex)
G4GLOB_DLL std::ostream G4cout
void StreamInfo(std::ostream &os) const

◆ EnablePolarisation()

◆ EnableSamplingTable()

G4bool G4EmParameters::EnableSamplingTable ( ) const

Definition at line 516 of file G4EmParameters.cc.

517{
518 return fSamplingTable;
519}

Referenced by G4SeltzerBergerModel::Initialise().

◆ FactorForAngleLimit()

G4double G4EmParameters::FactorForAngleLimit ( ) const

Definition at line 801 of file G4EmParameters.cc.

802{
803 return factorForAngleLimit;
804}

Referenced by G4WentzelOKandVIxSection::Initialise(), and G4CoulombScattering::InitialiseProcess().

◆ FillStepFunction()

void G4EmParameters::FillStepFunction ( const G4ParticleDefinition * part,
G4VEnergyLossProcess * proc ) const

Definition at line 1014 of file G4EmParameters.cc.

1015{
1016 fBParameters->FillStepFunction(part, proc);
1017}

◆ FluctuationType()

G4EmFluctuationType G4EmParameters::FluctuationType ( ) const

Definition at line 1093 of file G4EmParameters.cc.

1094{
1095 return fFluct;
1096}

Referenced by G4EmStandUtil::ModelOfFluctuations().

◆ Fluo()

G4bool G4EmParameters::Fluo ( ) const

Definition at line 253 of file G4EmParameters.cc.

254{
255 return fCParameters->Fluo();
256}

Referenced by G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), and G4VRadioactiveDecay::StreamInfo().

◆ FluoDirectory()

G4EmFluoDirectory G4EmParameters::FluoDirectory ( ) const

Definition at line 258 of file G4EmParameters.cc.

259{
260 return fCParameters->FluoDirectory();
261}

Referenced by G4AtomicTransitionManager::Initialise(), and StreamInfo().

◆ GeneralProcessActive()

G4bool G4EmParameters::GeneralProcessActive ( ) const

Definition at line 474 of file G4EmParameters.cc.

475{
476 return gener;
477}

Referenced by G4EmStandardPhysics_option4::ConstructProcess(), and G4EmStandardPhysicsSS::ConstructProcess().

◆ GetConversionType()

G4int G4EmParameters::GetConversionType ( ) const

Definition at line 1172 of file G4EmParameters.cc.

1173{
1174 return tripletConv;
1175}

Referenced by G4BetheHeitler5DModel::Initialise().

◆ GetDirectionalSplitting()

G4bool G4EmParameters::GetDirectionalSplitting ( ) const

Definition at line 1345 of file G4EmParameters.cc.

1345 {
1346 return fBParameters->GetDirectionalSplitting();
1347}

Referenced by G4EmBiasingManager::Initialise().

◆ GetDirectionalSplittingRadius()

G4double G4EmParameters::GetDirectionalSplittingRadius ( )

Definition at line 1372 of file G4EmParameters.cc.

1373{
1374 return fBParameters->GetDirectionalSplittingRadius();
1375}

Referenced by G4EmBiasingManager::Initialise().

◆ GetDirectionalSplittingTarget()

G4ThreeVector G4EmParameters::GetDirectionalSplittingTarget ( ) const

Definition at line 1361 of file G4EmParameters.cc.

1362{
1363 return fBParameters->GetDirectionalSplittingTarget();
1364}

Referenced by G4EmBiasingManager::Initialise().

◆ GetDirLEDATA()

◆ GetEmSaturation()

G4EmSaturation * G4EmParameters::GetEmSaturation ( )

Definition at line 592 of file G4EmParameters.cc.

593{
594 if (nullptr == emSaturation) {
595 G4AutoLock l(&emParametersMutex);
596 if (nullptr == emSaturation) {
597 emSaturation = new G4EmSaturation(1);
598 }
599 l.unlock();
600 }
601 birks = true;
602 return emSaturation;
603}
G4TemplateAutoLock< G4Mutex > G4AutoLock

◆ GetTimeStepModel()

G4ChemTimeStepModel G4EmParameters::GetTimeStepModel ( ) const

Definition at line 1623 of file G4EmParameters.cc.

1624{
1625 return fCParameters->GetChemTimeStepModel();
1626}

Referenced by G4EmDNAChemistry_option3::ConstructProcess(), G4EmDNAChemistry_option3::ConstructReactionTable(), and G4EmDNAChemistry_option3::ConstructTimeStepModel().

◆ Instance()

G4EmParameters * G4EmParameters::Instance ( )
static

Definition at line 71 of file G4EmParameters.cc.

72{
73 if(nullptr == theInstance) {
74 G4AutoLock l(&emParametersMutex);
75 if(nullptr == theInstance) {
76 static G4EmParameters manager;
77 theInstance = &manager;
78 }
79 l.unlock();
80 }
81 return theInstance;
82}
G4EmParameters(G4EmParameters &)=delete

Referenced by G4EmTableUtil::BuildMscProcess(), G4GammaConversionToMuons::BuildPhysicsTable(), G4TransportationWithMsc::BuildPhysicsTable(), G4EmBuilder::ConstructCharged(), G4EmBuilder::ConstructChargedSS(), G4EmDNABuilder::ConstructDNAProtonPhysics(), G4EmBuilder::ConstructElectronMscProcess(), G4EmBuilder::ConstructElectronSSProcess(), LBE::ConstructGeneral(), G4EmDNAChemistry_option3::ConstructProcess(), G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4RadioactiveDecayPhysics::ConstructProcess(), G4EmDNAChemistry_option3::ConstructReactionTable(), G4EmDNABuilder::ConstructStandardEmPhysics(), G4EmDNAChemistry_option3::ConstructTimeStepModel(), G4ECDecay::DecayIt(), G4ITDecay::DecayIt(), G4EmUtility::FindCrossSectionMax(), FTFP_BERT_TRV::FTFP_BERT_TRV(), G4CT_COUNT_IMPL(), G4DNABornIonisationModel::G4DNABornIonisationModel(), G4DNABornIonisationModel1::G4DNABornIonisationModel1(), G4DynamicParticleIonisation::G4DynamicParticleIonisation(), G4EmCalculator::G4EmCalculator(), G4EmDNAPhysics::G4EmDNAPhysics(), G4EmDNAPhysics_option1::G4EmDNAPhysics_option1(), G4EmDNAPhysics_option2::G4EmDNAPhysics_option2(), G4EmDNAPhysics_option3::G4EmDNAPhysics_option3(), G4EmDNAPhysics_option4::G4EmDNAPhysics_option4(), G4EmDNAPhysics_option5::G4EmDNAPhysics_option5(), G4EmDNAPhysics_option6::G4EmDNAPhysics_option6(), G4EmDNAPhysics_option7::G4EmDNAPhysics_option7(), G4EmDNAPhysics_option8::G4EmDNAPhysics_option8(), G4EmDNAPhysics_stationary::G4EmDNAPhysics_stationary(), G4EmDNAPhysics_stationary_option2::G4EmDNAPhysics_stationary_option2(), G4EmDNAPhysics_stationary_option4::G4EmDNAPhysics_stationary_option4(), G4EmDNAPhysics_stationary_option6::G4EmDNAPhysics_stationary_option6(), G4EmDNAPhysicsActivator::G4EmDNAPhysicsActivator(), G4EmElementXS::G4EmElementXS(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmModelActivator::G4EmModelActivator(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics::G4EmStandardPhysics(), G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), G4EmStandardPhysicsGS::G4EmStandardPhysicsGS(), G4EmStandardPhysicsSS::G4EmStandardPhysicsSS(), G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI(), G4LossTableBuilder::G4LossTableBuilder(), G4VEmAngularDistribution::G4VEmAngularDistribution(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4DNASolvationModelFactory::GetMacroDefinedModel(), G4AtomicTransitionManager::Initialise(), G4BetheBlochModel::Initialise(), G4BetheHeitler5DModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4BraggModel::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNABornIonisationModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAGeneralIonIonisationModel::Initialise(), G4DNARuddIonisationDynamicModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4EmBiasingManager::Initialise(), G4EmMultiModel::Initialise(), G4eplusTo2or3GammaModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4IonICRU73Data::Initialise(), G4IonParametrisedLossModel::Initialise(), G4MuPairProductionModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4RiGeMuPairProductionModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4UrbanMscModel::Initialise(), G4WentzelOKandVIxSection::Initialise(), G4WentzelVIModel::Initialise(), G4WentzelOKandVIxSection::InitialiseA(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), G4EmUtility::InitialiseElementSelectors(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4eIonisation::InitialiseEnergyLossProcess(), G4ePairProduction::InitialiseEnergyLossProcess(), G4hhIonisation::InitialiseEnergyLossProcess(), G4hIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4MuBremsstrahlung::InitialiseEnergyLossProcess(), G4MuIonisation::InitialiseEnergyLossProcess(), G4MuonToMuonPairProduction::InitialiseEnergyLossProcess(), G4MuPairProduction::InitialiseEnergyLossProcess(), G4PolarizedBremsstrahlung::InitialiseEnergyLossProcess(), G4PolarizedIonisation::InitialiseEnergyLossProcess(), G4UAtomicDeexcitation::InitialiseForNewRun(), G4PhysListUtil::InitialiseParameters(), G4VMscModel::InitialiseParameters(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4GammaGeneralProcess::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4DNASamplingTable::LoadData(), G4GoudsmitSaundersonTable::LoadMSCData(), G4CoulombScattering::MinPrimaryEnergy(), G4EmStandUtil::ModelOfFluctuations(), G4EmTableUtil::PrepareEmProcess(), G4EmTableUtil::PrepareMscProcess(), G4TransportationWithMsc::PreparePhysicsTable(), G4XrayReflection::ReadHenkeXrayData(), G4EmDataHandler::RetrievePhysicsTable(), G4MuPairProductionModel::RetrieveTables(), G4RiGeMuPairProductionModel::RetrieveTables(), G4VEmModel::SetLPMFlag(), G4VRadioactiveDecay::StreamInfo(), G4CoulombScattering::StreamProcessInfo(), and G4eBremsstrahlung::StreamProcessInfo().

◆ Integral()

G4bool G4EmParameters::Integral ( ) const

Definition at line 393 of file G4EmParameters.cc.

394{
395 return integral;
396}

◆ IsPrintLocked()

G4bool G4EmParameters::IsPrintLocked ( ) const

Definition at line 587 of file G4EmParameters.cc.

588{
589 return fIsPrinted;
590}

◆ LambdaFactor()

G4double G4EmParameters::LambdaFactor ( ) const

Definition at line 783 of file G4EmParameters.cc.

784{
785 return lambdaFactor;
786}

◆ LateralDisplacement()

G4bool G4EmParameters::LateralDisplacement ( ) const

Definition at line 338 of file G4EmParameters.cc.

339{
340 return lateralDisplacement;
341}

Referenced by G4VMscModel::InitialiseParameters().

◆ LateralDisplacementAlg96()

G4bool G4EmParameters::LateralDisplacementAlg96 ( ) const

Definition at line 349 of file G4EmParameters.cc.

350{
351 return lateralDisplacementAlg96;
352}

Referenced by G4UrbanMscModel::Initialise().

◆ LinearLossLimit()

G4double G4EmParameters::LinearLossLimit ( ) const

Definition at line 729 of file G4EmParameters.cc.

730{
731 return linLossLimit;
732}

◆ LivermoreDataDir()

const G4String & G4EmParameters::LivermoreDataDir ( )

Definition at line 1205 of file G4EmParameters.cc.

1206{
1207 return fCParameters->LivermoreDataDir();
1208}

◆ LossFluctuation()

G4bool G4EmParameters::LossFluctuation ( ) const

Definition at line 198 of file G4EmParameters.cc.

199{
200 return lossFluctuation;
201}

◆ LowestElectronEnergy()

G4double G4EmParameters::LowestElectronEnergy ( ) const

Definition at line 667 of file G4EmParameters.cc.

668{
669 return lowestElectronEnergy;
670}

◆ LowestMuHadEnergy()

G4double G4EmParameters::LowestMuHadEnergy ( ) const

Definition at line 678 of file G4EmParameters.cc.

679{
680 return lowestMuHadEnergy;
681}

◆ LowestTripletEnergy()

G4double G4EmParameters::LowestTripletEnergy ( ) const

Definition at line 689 of file G4EmParameters.cc.

690{
691 return lowestTripletEnergy;
692}

Referenced by G4eplusTo2or3GammaModel::Initialise().

◆ LPM()

G4bool G4EmParameters::LPM ( ) const

Definition at line 220 of file G4EmParameters.cc.

221{
222 return flagLPM;
223}

Referenced by G4eBremsstrahlungRelModel::Initialise(), and G4eBremsstrahlung::StreamProcessInfo().

◆ MaxDNAElectronEnergy()

G4double G4EmParameters::MaxDNAElectronEnergy ( ) const

Definition at line 974 of file G4EmParameters.cc.

975{
976 return fCParameters->MaxDNAElectronEnergy();
977}

Referenced by G4EmDNAPhysics::ConstructProcess().

◆ MaxDNAIonEnergy()

G4double G4EmParameters::MaxDNAIonEnergy ( ) const

Definition at line 985 of file G4EmParameters.cc.

986{
987 return fCParameters->MaxDNAIonEnergy();
988}

Referenced by G4EmDNAPhysics::ConstructProcess().

◆ MaxEnergyFor5DMuPair()

G4double G4EmParameters::MaxEnergyFor5DMuPair ( ) const

Definition at line 711 of file G4EmParameters.cc.

712{
713 return max5DEnergyForMuPair;
714}

Referenced by G4GammaConversionToMuons::BuildPhysicsTable().

◆ MaxEnergyForCSDARange()

G4double G4EmParameters::MaxEnergyForCSDARange ( ) const

Definition at line 656 of file G4EmParameters.cc.

657{
658 return maxKinEnergyCSDA;
659}

◆ MaxKinEnergy()

◆ MaxNIELEnergy()

◆ MinKinEnergy()

◆ MscEnergyLimit()

◆ MscGeomFactor()

G4double G4EmParameters::MscGeomFactor ( ) const

Definition at line 891 of file G4EmParameters.cc.

892{
893 return geomFactor;
894}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscLambdaLimit()

G4double G4EmParameters::MscLambdaLimit ( ) const

Definition at line 927 of file G4EmParameters.cc.

928{
929 return lambdaLimit;
930}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscMuHadRangeFactor()

G4double G4EmParameters::MscMuHadRangeFactor ( ) const

Definition at line 873 of file G4EmParameters.cc.

874{
875 return rangeFactorMuHad;
876}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscMuHadStepLimitType()

G4MscStepLimitType G4EmParameters::MscMuHadStepLimitType ( ) const

Definition at line 1126 of file G4EmParameters.cc.

1127{
1128 return mscStepLimitMuHad;
1129}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscPositronCorrection()

G4bool G4EmParameters::MscPositronCorrection ( ) const

Definition at line 532 of file G4EmParameters.cc.

533{
534 return fMscPosiCorr;
535}

Referenced by G4UrbanMscModel::Initialise().

◆ MscRangeFactor()

G4double G4EmParameters::MscRangeFactor ( ) const

Definition at line 855 of file G4EmParameters.cc.

856{
857 return rangeFactor;
858}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscSafetyFactor()

G4double G4EmParameters::MscSafetyFactor ( ) const

Definition at line 909 of file G4EmParameters.cc.

910{
911 return safetyFactor;
912}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscSkin()

G4double G4EmParameters::MscSkin ( ) const

Definition at line 945 of file G4EmParameters.cc.

946{
947 return skin;
948}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscStepLimitType()

G4MscStepLimitType G4EmParameters::MscStepLimitType ( ) const

Definition at line 1115 of file G4EmParameters.cc.

1116{
1117 return mscStepLimit;
1118}

Referenced by G4VMscModel::InitialiseParameters().

◆ MscThetaLimit()

G4double G4EmParameters::MscThetaLimit ( ) const

◆ MuHadBremsstrahlungTh()

G4double G4EmParameters::MuHadBremsstrahlungTh ( ) const

◆ MuHadLateralDisplacement()

G4bool G4EmParameters::MuHadLateralDisplacement ( ) const

Definition at line 360 of file G4EmParameters.cc.

361{
362 return muhadLateralDisplacement;
363}

Referenced by G4VMscModel::InitialiseParameters().

◆ NuclearFormfactorType()

G4NuclearFormfactorType G4EmParameters::NuclearFormfactorType ( ) const

Definition at line 1149 of file G4EmParameters.cc.

1150{
1151 return nucFormfactor;
1152}

Referenced by G4eSingleCoulombScatteringModel::Initialise(), and G4WentzelOKandVIxSection::Initialise().

◆ NumberForFreeVector()

G4int G4EmParameters::NumberForFreeVector ( ) const

Definition at line 1071 of file G4EmParameters.cc.

1072{
1073 return nForFreeVector;
1074}

◆ NumberOfBins()

G4int G4EmParameters::NumberOfBins ( ) const

Definition at line 1019 of file G4EmParameters.cc.

1020{
1021 return nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
1022}
int G4lrint(double ad)
Definition templates.hh:134

◆ NumberOfBinsPerDecade()

◆ OnIsolated()

G4bool G4EmParameters::OnIsolated ( ) const

Definition at line 505 of file G4EmParameters.cc.

506{
507 return onIsolated;
508}

Referenced by G4BetheHeitler5DModel::Initialise().

◆ operator=()

G4EmParameters & G4EmParameters::operator= ( const G4EmParameters & right)
delete

◆ ParticlesPAI()

const std::vector< G4String > & G4EmParameters::ParticlesPAI ( ) const

Definition at line 1223 of file G4EmParameters.cc.

1224{
1225 return fBParameters->ParticlesPAI();
1226}

◆ PhotoeffectBelowKShell()

G4bool G4EmParameters::PhotoeffectBelowKShell ( ) const

Definition at line 521 of file G4EmParameters.cc.

522{
523 return fPEKShell;
524}

Referenced by G4PEEffectFluoModel::Initialise().

◆ Pixe()

G4bool G4EmParameters::Pixe ( ) const

Definition at line 316 of file G4EmParameters.cc.

317{
318 return fCParameters->Pixe();
319}

Referenced by G4VAtomDeexcitation::InitialiseAtomicDeexcitation().

◆ PIXECrossSectionModel()

const G4String & G4EmParameters::PIXECrossSectionModel ( )

Definition at line 1183 of file G4EmParameters.cc.

1184{
1185 return fCParameters->PIXECrossSectionModel();
1186}

Referenced by G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), and G4UAtomicDeexcitation::InitialiseForNewRun().

◆ PIXEElectronCrossSectionModel()

const G4String & G4EmParameters::PIXEElectronCrossSectionModel ( )

Definition at line 1194 of file G4EmParameters.cc.

1195{
1196 return fCParameters->PIXEElectronCrossSectionModel();
1197}

Referenced by G4PenelopeIonisationModel::Initialise(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), and G4UAtomicDeexcitation::InitialiseForNewRun().

◆ PositronAtRestModelType()

G4PositronAtRestModelType G4EmParameters::PositronAtRestModelType ( ) const

Definition at line 1104 of file G4EmParameters.cc.

1105{
1106 return fPositronium;
1107}

◆ QuantumEntanglement()

G4bool G4EmParameters::QuantumEntanglement ( ) const

Definition at line 1334 of file G4EmParameters.cc.

1335{
1336 return fBParameters->QuantumEntanglement();
1337}

◆ RegionsDNA()

const std::vector< G4String > & G4EmParameters::RegionsDNA ( ) const

Definition at line 1256 of file G4EmParameters.cc.

1257{
1258 return fCParameters->RegionsDNA();
1259}

◆ RegionsMicroElec()

const std::vector< G4String > & G4EmParameters::RegionsMicroElec ( ) const

Definition at line 1244 of file G4EmParameters.cc.

1245{
1246 return fCParameters->RegionsMicroElec();
1247}

◆ RegionsPAI()

const std::vector< G4String > & G4EmParameters::RegionsPAI ( ) const

Definition at line 1228 of file G4EmParameters.cc.

1229{
1230 return fBParameters->RegionsPAI();
1231}

◆ RegionsPhysics()

const std::vector< G4String > & G4EmParameters::RegionsPhysics ( ) const

Definition at line 1272 of file G4EmParameters.cc.

1273{
1274 return fBParameters->RegionsPhysics();
1275}

◆ RetrieveMuDataFromFile()

G4bool G4EmParameters::RetrieveMuDataFromFile ( ) const

Definition at line 489 of file G4EmParameters.cc.

490{
491 return fMuDataFromFile;
492}

Referenced by G4MuPairProductionModel::Initialise(), and G4RiGeMuPairProductionModel::Initialise().

◆ ScreeningFactor()

G4double G4EmParameters::ScreeningFactor ( ) const

Definition at line 963 of file G4EmParameters.cc.

964{
965 return factorScreen;
966}

Referenced by G4WentzelOKandVIxSection::InitialiseA().

◆ Set3GammaAnnihilationOnFly()

void G4EmParameters::Set3GammaAnnihilationOnFly ( G4bool v)

Definition at line 559 of file G4EmParameters.cc.

560{
561 if(IsLocked()) { return; }
562 f3GammaAnnihilationOnFly = v;
563}

Referenced by G4EmStandardPhysics_option4::G4EmStandardPhysics_option4().

◆ SetANSTOFluoDir()

void G4EmParameters::SetANSTOFluoDir ( G4bool val)

Definition at line 275 of file G4EmParameters.cc.

276{
277 if(IsLocked()) { return; }
278 fCParameters->SetANSTOFluoDir(val);
279}

◆ SetApplyCuts()

void G4EmParameters::SetApplyCuts ( G4bool val)

Definition at line 236 of file G4EmParameters.cc.

237{
238 if(IsLocked()) { return; }
239 applyCuts = val;
240}

Referenced by G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), and G4EmStandardPhysics_option2::G4EmStandardPhysics_option2().

◆ SetAuger()

void G4EmParameters::SetAuger ( G4bool val)

Definition at line 287 of file G4EmParameters.cc.

288{
289 if(IsLocked()) { return; }
290 fCParameters->SetAuger(val);
291}

Referenced by G4RadioactiveDecayPhysics::ConstructProcess(), G4EmDNAPhysics::G4EmDNAPhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmStandardPhysicsSS::G4EmStandardPhysicsSS(), and SetAugerCascade().

◆ SetAugerCascade()

void G4EmParameters::SetAugerCascade ( G4bool val)
inline

Definition at line 148 of file G4EmParameters.hh.

148{ SetAuger(val); };
void SetAuger(G4bool val)

Referenced by LBE::ConstructGeneral().

◆ SetBeardenFluoDir()

void G4EmParameters::SetBeardenFluoDir ( G4bool val)

Definition at line 269 of file G4EmParameters.cc.

270{
271 if(IsLocked()) { return; }
272 fCParameters->SetBeardenFluoDir(val);
273}

◆ SetBirksActive()

void G4EmParameters::SetBirksActive ( G4bool val)

Definition at line 409 of file G4EmParameters.cc.

410{
411 if(IsLocked()) { return; }
412 birks = val;
413 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
414}

◆ SetBremsstrahlungTh()

void G4EmParameters::SetBremsstrahlungTh ( G4double val)

Definition at line 734 of file G4EmParameters.cc.

735{
736 if(IsLocked()) { return; }
737 if(val > 0.0) {
738 bremsTh = val;
739 } else {
741 ed << "Value of bremsstrahlung threshold is out of range: "
742 << val/GeV << " GeV is ignored";
743 PrintWarning(ed);
744 }
745}
std::ostringstream G4ExceptionDescription

◆ SetBuildCSDARange()

void G4EmParameters::SetBuildCSDARange ( G4bool val)

Definition at line 203 of file G4EmParameters.cc.

204{
205 if(IsLocked()) { return; }
206 buildCSDARange = val;
207}

◆ SetConversionType()

void G4EmParameters::SetConversionType ( G4int val)

Definition at line 1166 of file G4EmParameters.cc.

1167{
1168 if(IsLocked()) { return; }
1169 tripletConv = val;
1170}

◆ SetDeexActiveRegion()

void G4EmParameters::SetDeexActiveRegion ( const G4String & region,
G4bool fdeex,
G4bool fauger,
G4bool fpixe )

Definition at line 1289 of file G4EmParameters.cc.

1291{
1292 if(IsLocked()) { return; }
1293 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1294}

◆ SetDeexcitationIgnoreCut()

void G4EmParameters::SetDeexcitationIgnoreCut ( G4bool val)

Definition at line 321 of file G4EmParameters.cc.

322{
323 if(IsLocked()) { return; }
324 fCParameters->SetDeexcitationIgnoreCut(val);
325}

Referenced by G4RadioactiveDecayPhysics::ConstructProcess(), and G4EmDNAPhysics::G4EmDNAPhysics().

◆ SetDefaults()

◆ SetDirectionalSplitting()

void G4EmParameters::SetDirectionalSplitting ( G4bool v)

Definition at line 1349 of file G4EmParameters.cc.

1350{
1351 if(IsLocked()) { return; }
1352 fBParameters->SetDirectionalSplitting(v);
1353}

◆ SetDirectionalSplittingRadius()

void G4EmParameters::SetDirectionalSplittingRadius ( G4double r)

Definition at line 1366 of file G4EmParameters.cc.

1367{
1368 if(IsLocked()) { return; }
1369 fBParameters->SetDirectionalSplittingRadius(r);
1370}

◆ SetDirectionalSplittingTarget()

void G4EmParameters::SetDirectionalSplittingTarget ( const G4ThreeVector & v)

Definition at line 1355 of file G4EmParameters.cc.

1356{
1357 if(IsLocked()) { return; }
1358 fBParameters->SetDirectionalSplittingTarget(v);
1359}

◆ SetDNAElectronMsc()

void G4EmParameters::SetDNAElectronMsc ( G4bool val)

Definition at line 456 of file G4EmParameters.cc.

457{
458 if(IsLocked()) { return; }
459 fCParameters->SetDNAElectronMsc(val);
460 if(val) { ActivateDNA(); }
461}

◆ SetDNAeSolvationSubType()

void G4EmParameters::SetDNAeSolvationSubType ( G4DNAModelSubType val)

Definition at line 1154 of file G4EmParameters.cc.

1155{
1156 if(IsLocked()) { return; }
1157 fCParameters->SetDNAeSolvationSubType(val);
1158 ActivateDNA();
1159}

◆ SetDNAFast()

◆ SetDNAStationary()

◆ SetEmSaturation()

void G4EmParameters::SetEmSaturation ( G4EmSaturation * ptr)

Definition at line 479 of file G4EmParameters.cc.

480{
481 if(IsLocked()) { return; }
482 birks = (nullptr != ptr);
483 if(emSaturation != ptr) {
484 delete emSaturation;
485 emSaturation = ptr;
486 }
487}

◆ SetEnablePolarisation()

void G4EmParameters::SetEnablePolarisation ( G4bool val)

Definition at line 398 of file G4EmParameters.cc.

399{
400 if(IsLocked()) { return; }
401 fPolarisation = val;
402}

Referenced by G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics().

◆ SetEnableSamplingTable()

void G4EmParameters::SetEnableSamplingTable ( G4bool val)

Definition at line 510 of file G4EmParameters.cc.

511{
512 if(IsLocked()) { return; }
513 fSamplingTable = val;
514}

◆ SetFactorForAngleLimit()

void G4EmParameters::SetFactorForAngleLimit ( G4double val)

Definition at line 788 of file G4EmParameters.cc.

789{
790 if(IsLocked()) { return; }
791 if(val > 0.0) {
792 factorForAngleLimit = val;
793 } else {
795 ed << "Value of factor for enegry limit is out of range: "
796 << val << " is ignored";
797 PrintWarning(ed);
798 }
799}

◆ SetFluctuationsForRegion()

void G4EmParameters::SetFluctuationsForRegion ( const G4String & regionName,
G4bool flag )

Definition at line 1387 of file G4EmParameters.cc.

1388{
1389 if (IsLocked()) { return; }
1390 G4String ss = fBParameters->CheckRegion(nam);
1391 if (!fluctRegions.empty()) {
1392 for (auto const& p : fluctRegions) { if (p.first == ss) { return; } }
1393 }
1394 fluctRegions.push_back(std::make_pair(ss, flag));
1395}

◆ SetFluctuationType()

◆ SetFluo()

◆ SetFluoDirectory()

void G4EmParameters::SetFluoDirectory ( G4EmFluoDirectory val)

Definition at line 263 of file G4EmParameters.cc.

264{
265 if(IsLocked()) { return; }
266 fCParameters->SetFluoDirectory(val);
267}

◆ SetGeneralProcessActive()

void G4EmParameters::SetGeneralProcessActive ( G4bool val)

◆ SetIntegral()

void G4EmParameters::SetIntegral ( G4bool val)

Definition at line 387 of file G4EmParameters.cc.

388{
389 if(IsLocked()) { return; }
390 integral = val;
391}

◆ SetIsPrintedFlag()

void G4EmParameters::SetIsPrintedFlag ( G4bool val)

Definition at line 582 of file G4EmParameters.cc.

583{
584 fIsPrinted = val;
585}

◆ SetLambdaFactor()

void G4EmParameters::SetLambdaFactor ( G4double val)

Definition at line 770 of file G4EmParameters.cc.

771{
772 if(IsLocked()) { return; }
773 if(val > 0.0 && val < 1.0) {
774 lambdaFactor = val;
775 } else {
777 ed << "Value of lambda factor is out of range: " << val
778 << " is ignored";
779 PrintWarning(ed);
780 }
781}

◆ SetLateralDisplacement()

void G4EmParameters::SetLateralDisplacement ( G4bool val)

Definition at line 332 of file G4EmParameters.cc.

333{
334 if(IsLocked()) { return; }
335 lateralDisplacement = val;
336}

Referenced by G4EmStandardPhysics_option2::G4EmStandardPhysics_option2().

◆ SetLateralDisplacementAlg96()

void G4EmParameters::SetLateralDisplacementAlg96 ( G4bool val)

Definition at line 343 of file G4EmParameters.cc.

344{
345 if(IsLocked()) { return; }
346 lateralDisplacementAlg96 = val;
347}

Referenced by G4EmStandardPhysics_option3::G4EmStandardPhysics_option3().

◆ SetLinearLossLimit()

void G4EmParameters::SetLinearLossLimit ( G4double val)

Definition at line 716 of file G4EmParameters.cc.

717{
718 if(IsLocked()) { return; }
719 if(val > 0.0 && val < 1.0) {
720 linLossLimit = val;
721 } else {
723 ed << "Value of linLossLimit is out of range: " << val
724 << " is ignored";
725 PrintWarning(ed);
726 }
727}

◆ SetLivermoreDataDir()

void G4EmParameters::SetLivermoreDataDir ( const G4String & sss)

Definition at line 1199 of file G4EmParameters.cc.

1200{
1201 if(IsLocked()) { return; }
1202 fCParameters->SetLivermoreDataDir(sss);
1203}

◆ SetLossFluctuations()

void G4EmParameters::SetLossFluctuations ( G4bool val)

Definition at line 192 of file G4EmParameters.cc.

193{
194 if(IsLocked()) { return; }
195 lossFluctuation = val;
196}

◆ SetLowestElectronEnergy()

◆ SetLowestMuHadEnergy()

void G4EmParameters::SetLowestMuHadEnergy ( G4double val)

Definition at line 672 of file G4EmParameters.cc.

673{
674 if(IsLocked()) { return; }
675 if(val >= 0.0) { lowestMuHadEnergy = val; }
676}

◆ SetLowestTripletEnergy()

void G4EmParameters::SetLowestTripletEnergy ( G4double val)

Definition at line 683 of file G4EmParameters.cc.

684{
685 if(IsLocked()) { return; }
686 if(val > 0.0) { lowestTripletEnergy = val; }
687}

◆ SetLPM()

void G4EmParameters::SetLPM ( G4bool val)

Definition at line 214 of file G4EmParameters.cc.

215{
216 if(IsLocked()) { return; }
217 flagLPM = val;
218}

◆ SetMaxDNAElectronEnergy()

void G4EmParameters::SetMaxDNAElectronEnergy ( G4double val)

Definition at line 968 of file G4EmParameters.cc.

969{
970 if (IsLocked()) { return; }
971 fCParameters->SetMaxDNAElectronEnergy(val);
972}

◆ SetMaxDNAIonEnergy()

void G4EmParameters::SetMaxDNAIonEnergy ( G4double val)

Definition at line 979 of file G4EmParameters.cc.

980{
981 if (IsLocked()) { return; }
982 fCParameters->SetMaxDNAIonEnergy(val);
983}

◆ SetMaxEnergy()

void G4EmParameters::SetMaxEnergy ( G4double val)

Definition at line 623 of file G4EmParameters.cc.

624{
625 if(IsLocked()) { return; }
626 if(val > std::max(minKinEnergy,599.9*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
627 maxKinEnergy = val;
628 } else {
630 ed << "Value of MaxKinEnergy is out of range: "
631 << val/CLHEP::GeV
632 << " GeV is ignored; allowed range 600 MeV - 1.e+7 TeV";
633 PrintWarning(ed);
634 }
635}

Referenced by G4EmDNAPhysics::G4EmDNAPhysics().

◆ SetMaxEnergyFor5DMuPair()

void G4EmParameters::SetMaxEnergyFor5DMuPair ( G4double val)

Definition at line 705 of file G4EmParameters.cc.

706{
707 if(IsLocked()) { return; }
708 if(val > 0.0) { max5DEnergyForMuPair = val; }
709}

◆ SetMaxEnergyForCSDARange()

void G4EmParameters::SetMaxEnergyForCSDARange ( G4double val)

Definition at line 642 of file G4EmParameters.cc.

643{
644 if(IsLocked()) { return; }
645 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
646 maxKinEnergyCSDA = val;
647 } else {
649 ed << "Value of MaxKinEnergyCSDA is out of range: "
650 << val/CLHEP::GeV << " GeV is ignored; allowed range "
651 << minKinEnergy << " MeV - 100 TeV";
652 PrintWarning(ed);
653 }
654}

◆ SetMaxNIELEnergy()

void G4EmParameters::SetMaxNIELEnergy ( G4double val)

◆ SetMinEnergy()

void G4EmParameters::SetMinEnergy ( G4double val)

Definition at line 605 of file G4EmParameters.cc.

606{
607 if(IsLocked()) { return; }
608 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
609 minKinEnergy = val;
610 } else {
612 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
613 << " MeV is ignored";
614 PrintWarning(ed);
615 }
616}

Referenced by G4EmDNAPhysics::G4EmDNAPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), and G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI().

◆ SetMscEnergyLimit()

void G4EmParameters::SetMscEnergyLimit ( G4double val)

Definition at line 824 of file G4EmParameters.cc.

825{
826 if(IsLocked()) { return; }
827 if(val >= 0.0) {
828 energyLimit = val;
829 } else {
831 ed << "Value of msc energy limit is out of range: "
832 << val << " is ignored";
833 PrintWarning(ed);
834 }
835}

Referenced by G4EmStandardPhysicsGS::G4EmStandardPhysicsGS().

◆ SetMscGeomFactor()

void G4EmParameters::SetMscGeomFactor ( G4double val)

Definition at line 878 of file G4EmParameters.cc.

879{
880 if(IsLocked()) { return; }
881 if(val >= 1.0) {
882 geomFactor = val;
883 } else {
885 ed << "Value of geomFactor is out of range: "
886 << val << " is ignored";
887 PrintWarning(ed);
888 }
889}

◆ SetMscLambdaLimit()

void G4EmParameters::SetMscLambdaLimit ( G4double val)

Definition at line 914 of file G4EmParameters.cc.

915{
916 if(IsLocked()) { return; }
917 if(val >= 0.0) {
918 lambdaLimit = val;
919 } else {
921 ed << "Value of lambdaLimit is out of range: "
922 << val << " is ignored";
923 PrintWarning(ed);
924 }
925}

◆ SetMscMuHadRangeFactor()

void G4EmParameters::SetMscMuHadRangeFactor ( G4double val)

Definition at line 860 of file G4EmParameters.cc.

861{
862 if(IsLocked()) { return; }
863 if(val > 0.0 && val < 1.0) {
864 rangeFactorMuHad = val;
865 } else {
867 ed << "Value of rangeFactorMuHad is out of range: "
868 << val << " is ignored";
869 PrintWarning(ed);
870 }
871}

◆ SetMscMuHadStepLimitType()

void G4EmParameters::SetMscMuHadStepLimitType ( G4MscStepLimitType val)

Definition at line 1120 of file G4EmParameters.cc.

1121{
1122 if(IsLocked()) { return; }
1123 mscStepLimitMuHad = val;
1124}

◆ SetMscPositronCorrection()

void G4EmParameters::SetMscPositronCorrection ( G4bool v)

Definition at line 537 of file G4EmParameters.cc.

538{
539 if(IsLocked()) { return; }
540 fMscPosiCorr = v;
541}

◆ SetMscRangeFactor()

void G4EmParameters::SetMscRangeFactor ( G4double val)

Definition at line 842 of file G4EmParameters.cc.

843{
844 if(IsLocked()) { return; }
845 if(val > 0.0 && val < 1.0) {
846 rangeFactor = val;
847 } else {
849 ed << "Value of rangeFactor is out of range: "
850 << val << " is ignored";
851 PrintWarning(ed);
852 }
853}

Referenced by G4EmDNAPhysics::G4EmDNAPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), and G4EmStandardPhysicsGS::G4EmStandardPhysicsGS().

◆ SetMscSafetyFactor()

void G4EmParameters::SetMscSafetyFactor ( G4double val)

Definition at line 896 of file G4EmParameters.cc.

897{
898 if(IsLocked()) { return; }
899 if(val >= 0.1) {
900 safetyFactor = val;
901 } else {
903 ed << "Value of safetyFactor is out of range: "
904 << val << " is ignored";
905 PrintWarning(ed);
906 }
907}

◆ SetMscSkin()

void G4EmParameters::SetMscSkin ( G4double val)

Definition at line 932 of file G4EmParameters.cc.

933{
934 if(IsLocked()) { return; }
935 if(val >= 1.0) {
936 skin = val;
937 } else {
939 ed << "Value of skin is out of range: "
940 << val << " is ignored";
941 PrintWarning(ed);
942 }
943}

Referenced by G4EmDNAPhysics::G4EmDNAPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), and G4EmStandardPhysicsGS::G4EmStandardPhysicsGS().

◆ SetMscStepLimitType()

◆ SetMscThetaLimit()

void G4EmParameters::SetMscThetaLimit ( G4double val)

Definition at line 806 of file G4EmParameters.cc.

807{
808 if(IsLocked()) { return; }
809 if(val >= 0.0 && val <= pi) {
810 thetaLimit = val;
811 } else {
813 ed << "Value of polar angle limit is out of range: "
814 << val << " is ignored";
815 PrintWarning(ed);
816 }
817}

Referenced by G4EmStandardPhysicsSS::G4EmStandardPhysicsSS(), and G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI().

◆ SetMuHadBremsstrahlungTh()

void G4EmParameters::SetMuHadBremsstrahlungTh ( G4double val)

Definition at line 752 of file G4EmParameters.cc.

753{
754 if(IsLocked()) { return; }
755 if(val > 0.0) {
756 bremsMuHadTh = val;
757 } else {
759 ed << "Value of bremsstrahlung threshold is out of range: "
760 << val/GeV << " GeV is ignored";
761 PrintWarning(ed);
762 }
763}

◆ SetMuHadLateralDisplacement()

◆ SetNuclearFormfactorType()

void G4EmParameters::SetNuclearFormfactorType ( G4NuclearFormfactorType val)

Definition at line 1143 of file G4EmParameters.cc.

1144{
1145 if(IsLocked()) { return; }
1146 nucFormfactor = val;
1147}

◆ SetNumberForFreeVector()

void G4EmParameters::SetNumberForFreeVector ( G4int val)

Definition at line 1065 of file G4EmParameters.cc.

1066{
1067 if(IsLocked()) { return; }
1068 nForFreeVector = val;
1069}

◆ SetNumberOfBinsPerDecade()

void G4EmParameters::SetNumberOfBinsPerDecade ( G4int val)

Definition at line 1024 of file G4EmParameters.cc.

1025{
1026 if(IsLocked()) { return; }
1027 if(val >= 5 && val < 1000000) {
1028 nbinsPerDecade = val;
1029 } else {
1031 ed << "Value of number of bins per decade is out of range: "
1032 << val << " is ignored";
1033 PrintWarning(ed);
1034 }
1035}

Referenced by G4EmDNAPhysics::G4EmDNAPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), G4EmStandardPhysicsGS::G4EmStandardPhysicsGS(), and G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI().

◆ SetOnIsolated()

void G4EmParameters::SetOnIsolated ( G4bool val)

Definition at line 499 of file G4EmParameters.cc.

500{
501 if(IsLocked()) { return; }
502 onIsolated = val;
503}

◆ SetPhotoeffectBelowKShell()

void G4EmParameters::SetPhotoeffectBelowKShell ( G4bool v)

Definition at line 526 of file G4EmParameters.cc.

527{
528 if(IsLocked()) { return; }
529 fPEKShell = v;
530}

◆ SetPixe()

void G4EmParameters::SetPixe ( G4bool val)

Definition at line 310 of file G4EmParameters.cc.

311{
312 if(IsLocked()) { return; }
313 fCParameters->SetPixe(val);
314}

Referenced by G4EmStandardPhysicsSS::G4EmStandardPhysicsSS().

◆ SetPIXECrossSectionModel()

void G4EmParameters::SetPIXECrossSectionModel ( const G4String & sss)

Definition at line 1177 of file G4EmParameters.cc.

1178{
1179 if(IsLocked()) { return; }
1180 fCParameters->SetPIXECrossSectionModel(sss);
1181}

◆ SetPIXEElectronCrossSectionModel()

void G4EmParameters::SetPIXEElectronCrossSectionModel ( const G4String & sss)

Definition at line 1188 of file G4EmParameters.cc.

1189{
1190 if(IsLocked()) { return; }
1191 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1192}

Referenced by G4EmPenelopePhysics::G4EmPenelopePhysics().

◆ SetPositronAtRestModelType()

◆ SetProcessBiasingFactor()

void G4EmParameters::SetProcessBiasingFactor ( const G4String & procname,
G4double val,
G4bool wflag )

Definition at line 1297 of file G4EmParameters.cc.

1299{
1300 if(IsLocked()) { return; }
1301 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1302}

◆ SetQuantumEntanglement()

void G4EmParameters::SetQuantumEntanglement ( G4bool v)

Definition at line 1339 of file G4EmParameters.cc.

1340{
1341 if(IsLocked()) { return; }
1342 fBParameters->SetQuantumEntanglement(v);
1343}

◆ SetRetrieveMuDataFromFile()

void G4EmParameters::SetRetrieveMuDataFromFile ( G4bool v)

Definition at line 494 of file G4EmParameters.cc.

495{
496 fMuDataFromFile = v;
497}

◆ SetScreeningFactor()

void G4EmParameters::SetScreeningFactor ( G4double val)

Definition at line 950 of file G4EmParameters.cc.

951{
952 if(IsLocked()) { return; }
953 if(val > 0.0) {
954 factorScreen = val;
955 } else {
957 ed << "Value of factorScreen is out of range: "
958 << val << " is ignored";
959 PrintWarning(ed);
960 }
961}

◆ SetSingleScatteringType()

void G4EmParameters::SetSingleScatteringType ( G4eSingleScatteringType val)

Definition at line 1131 of file G4EmParameters.cc.

1132{
1133 if(IsLocked()) { return; }
1134 fSStype = val;
1135}

◆ SetStepFunction()

◆ SetStepFunctionIons()

◆ SetStepFunctionLightIons()

◆ SetStepFunctionMuHad()

◆ SetSubCutRegion()

void G4EmParameters::SetSubCutRegion ( const G4String & region = "")

Definition at line 1282 of file G4EmParameters.cc.

1283{
1284 if(IsLocked()) { return; }
1285 fBParameters->SetSubCutRegion(region);
1286}

◆ SetTimeStepModel()

void G4EmParameters::SetTimeStepModel ( const G4ChemTimeStepModel & model)

Definition at line 1618 of file G4EmParameters.cc.

1619{
1620 fCParameters-> SetChemTimeStepModel(model);
1621}

◆ SetTransportationWithMsc()

void G4EmParameters::SetTransportationWithMsc ( G4TransportationWithMscType val)

Definition at line 1076 of file G4EmParameters.cc.

1077{
1078 if(IsLocked()) { return; }
1079 fTransportationWithMsc = val;
1080}

Referenced by G4EmStandardPhysics_option1::G4EmStandardPhysics_option1().

◆ SetUseCutAsFinalRange()

void G4EmParameters::SetUseCutAsFinalRange ( G4bool val)

Definition at line 225 of file G4EmParameters.cc.

226{
227 if(IsLocked()) { return; }
228 cutAsFinalRange = val;
229}

◆ SetUseEPICS2017XS()

void G4EmParameters::SetUseEPICS2017XS ( G4bool v)

Definition at line 548 of file G4EmParameters.cc.

549{
550 if(IsLocked()) { return; }
551 fUseEPICS2017XS = v;
552}

◆ SetUseICRU90Data()

◆ SetUseMottCorrection()

◆ SetUseRiGePairProductionModel()

void G4EmParameters::SetUseRiGePairProductionModel ( G4bool v)

Definition at line 570 of file G4EmParameters.cc.

571{
572 if (IsLocked()) { return; }
573 fUseRiGePairProductionModel = v;
574}

◆ SetVerbose()

◆ SetWorkerVerbose()

void G4EmParameters::SetWorkerVerbose ( G4int val)

Definition at line 1054 of file G4EmParameters.cc.

1055{
1056 if(IsLocked()) { return; }
1057 workerVerbose = val;
1058}

◆ SetXDB_EADLFluoDir()

void G4EmParameters::SetXDB_EADLFluoDir ( G4bool val)

Definition at line 281 of file G4EmParameters.cc.

282{
283 if(IsLocked()) { return; }
284 fCParameters->SetXDB_EADLFluoDir(val);
285}

◆ SingleScatteringType()

G4eSingleScatteringType G4EmParameters::SingleScatteringType ( ) const

Definition at line 1137 of file G4EmParameters.cc.

1138{
1139 return fSStype;
1140}

◆ StreamInfo()

void G4EmParameters::StreamInfo ( std::ostream & os) const

Definition at line 1403 of file G4EmParameters.cc.

1404{
1405 G4long prec = os.precision(5);
1406 os << "=======================================================================" << "\n";
1407 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1408 os << "=======================================================================" << "\n";
1409 os << "LPM effect enabled " <<flagLPM << "\n";
1410 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1411 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1412 const char* transportationWithMsc = "Disabled";
1413 if(fTransportationWithMsc == G4TransportationWithMscType::fEnabled) {
1414 transportationWithMsc = "Enabled";
1415 } else if (fTransportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
1416 transportationWithMsc = "MultipleSteps";
1417 }
1418 os << "Use combined TransportationWithMsc " <<transportationWithMsc << "\n";
1419 os << "Use general process " <<gener << "\n";
1420 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1421 os << "Enable photoeffect sampling below K-shell " <<fPEKShell << "\n";
1422 os << "Enable sampling of quantum entanglement "
1423 <<fBParameters->QuantumEntanglement() << "\n";
1424 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1425 os << "Min kinetic energy for tables "
1426 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1427 os << "Max kinetic energy for tables "
1428 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1429 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1430 os << "Verbose level " <<verbose << "\n";
1431 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1432 os << "Bremsstrahlung energy threshold above which \n"
1433 << " primary e+- is added to the list of secondary "
1434 <<G4BestUnit(bremsTh,"Energy") << "\n";
1435 os << "Bremsstrahlung energy threshold above which primary\n"
1436 << " muon/hadron is added to the list of secondary "
1437 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1438 G4String name3g = "SimplePositronium";
1439 if (fPositronium == fAllisonPositronium) { name3g = "AllisonPositronium"; }
1440 else if (fPositronium == fOrePowell) { name3g = "OrePowell"; }
1441 else if (fPositronium == fOrePowellPolar) { name3g = "OrePowellPolar"; }
1442 os << "Positron annihilation at rest model " << name3g << "\n";
1443
1444 os << "Enable 3 gamma annihilation on fly "
1445 << f3GammaAnnihilationOnFly << "\n";
1446 os << "Lowest triplet kinetic energy "
1447 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1448 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1449 os << "5D gamma conversion model type " <<tripletConv << "\n";
1450 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1451 if(max5DEnergyForMuPair>0.0) {
1452 os << "5D gamma conversion limit for muon pair "
1453 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1454 }
1455 os << "Use RiGe 5D e+e- pair production model by muons "
1456 << fUseRiGePairProductionModel << "\n";
1457 os << "Livermore data directory "
1458 << fCParameters->LivermoreDataDir() << "\n";
1459
1460 os << "=======================================================================" << "\n";
1461 os << "====== Ionisation Parameters ========" << "\n";
1462 os << "=======================================================================" << "\n";
1463 os << "Step function for e+- "
1464 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1465 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1466 os << "Step function for muons/hadrons "
1467 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1468 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1469 os << "Step function for light ions "
1470 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1471 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1472 os << "Step function for general ions "
1473 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1474 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1475 os << "Lowest e+e- kinetic energy "
1476 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1477 os << "Lowest muon/hadron kinetic energy "
1478 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1479 os << "Use ICRU90 data " << fICRU90 << "\n";
1480 os << "Fluctuations of dE/dx are enabled " << lossFluctuation << "\n";
1481 G4String namef = "Universal";
1482 if(fFluct == fUrbanFluctuation) { namef = "Urban"; }
1483 else if(fFluct == fDummyFluctuation) { namef = "Dummy"; }
1484 os << "Type of fluctuation model for leptons and hadrons " << namef << "\n";
1485 if (!fluctRegions.empty()) {
1486 if (lossFluctuation) {
1487 os << "Fluctuations of dE/dx are disabled in G4Regions: ";
1488 } else {
1489 os << "Fluctuations of dE/dx are enabled in G4Regions: ";
1490 }
1491 G4int n = 0;
1492 for (auto const& p : fluctRegions) {
1493 if (p.second != lossFluctuation) {
1494 os << p.first << " ";
1495 if (n > 0 && (n/2)*2 == n) {
1496 os << "\n" << " ";
1497 }
1498 ++n;
1499 }
1500 }
1501 os << "\n";
1502 }
1503 os << "Use built-in Birks saturation " << birks << "\n";
1504 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1505 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1506 os << "Enable angular generator interface "
1507 <<useAngGeneratorForIonisation << "\n";
1508 os << "Max kinetic energy for CSDA tables "
1509 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1510 os << "Max kinetic energy for NIEL computation "
1511 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1512 os << "Linear loss limit " <<linLossLimit << "\n";
1513 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1514
1515 os << "=======================================================================" << "\n";
1516 os << "====== Multiple Scattering Parameters ========" << "\n";
1517 os << "=======================================================================" << "\n";
1518 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1519 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1520 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1521 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1522 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1523 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1524 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1525 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1526 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1527 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1528 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1529 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1530 os << "Factor used for dynamic computation of angular \n"
1531 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1532 os << "Fixed angular limit between single \n"
1533 << " and multiple scattering "
1534 << thetaLimit/CLHEP::rad << " rad\n";
1535 os << "Upper energy limit for e+- multiple scattering "
1536 << energyLimit/CLHEP::MeV << " MeV\n";
1537 os << "Type of electron single scattering model " <<fSStype << "\n";
1538 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1539 os << "Screening factor " <<factorScreen << "\n";
1540 os << "=======================================================================" << "\n";
1541
1542 if(fCParameters->Fluo()) {
1543 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1544 os << "=======================================================================" << "\n";
1545 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1546 G4String named = "fluor";
1548 if(fdir == fluoBearden) { named = "fluor_Bearden"; }
1549 else if(fdir == fluoANSTO) { named = "fluor_ANSTO"; }
1550 else if(fdir == fluoXDB_EADL) { named = "fluor_XDB_EADL"; }
1551 os << "Directory in G4LEDATA for fluorescence data files " << named << "\n";
1552 os << "Auger electron cascade enabled "
1553 <<fCParameters->Auger() << "\n";
1554 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1555 os << "De-excitation module ignores cuts "
1556 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1557 os << "Type of PIXE cross section for hadrons "
1558 <<fCParameters->PIXECrossSectionModel() << "\n";
1559 os << "Type of PIXE cross section for e+- "
1560 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1561 os << "=======================================================================" << "\n";
1562 }
1563 if(fDNA) {
1564 os << "====== DNA Physics Parameters ========" << "\n";
1565 os << "=======================================================================" << "\n";
1566 os << "Max DNA electron energy "
1567 << fCParameters->MaxDNAElectronEnergy()/CLHEP::MeV << " MeV\n";
1568 os << "Max DNA ion energy "
1569 << fCParameters->MaxDNAIonEnergy()/CLHEP::MeV << " MeV\n";
1570 os << "Use fast sampling in DNA models "
1571 << fCParameters->DNAFast() << "\n";
1572 os << "Use Stationary option in DNA models "
1573 << fCParameters->DNAStationary() << "\n";
1574 os << "Use DNA with multiple scattering of e- "
1575 << fCParameters->DNAElectronMsc() << "\n";
1576 os << "Use DNA e- solvation model type "
1577 << fCParameters->DNAeSolvationSubType() << "\n";
1578 auto chemModel = fCParameters->GetChemTimeStepModel();
1579 if(fCParameters->GetChemTimeStepModel() != G4ChemTimeStepModel::Unknown)
1580 {
1581 std::vector<G4String> ChemModel{"Unknown","SBS","IRT","IRT_syn"};
1582 os << "Use DNA Chemistry model "
1583 << ChemModel.at((std::size_t)chemModel) << "\n";
1584 }
1585 os << "=======================================================================" << G4endl;
1586 }
1587 os.precision(prec);
1588}
G4EmFluoDirectory
@ fluoXDB_EADL
@ fUrbanFluctuation
@ fDummyFluctuation
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
#define G4BestUnit(a, b)
long G4long
Definition G4Types.hh:87
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4EmFluoDirectory FluoDirectory() const

Referenced by Dump(), and operator<<.

◆ TransportationWithMsc()

G4TransportationWithMscType G4EmParameters::TransportationWithMsc ( ) const

Definition at line 1082 of file G4EmParameters.cc.

1083{
1084 return fTransportationWithMsc;
1085}

Referenced by G4EmBuilder::ConstructElectronMscProcess(), and G4EmBuilder::ConstructElectronSSProcess().

◆ TypesDNA()

const std::vector< G4String > & G4EmParameters::TypesDNA ( ) const

Definition at line 1261 of file G4EmParameters.cc.

1262{
1263 return fCParameters->TypesDNA();
1264}

◆ TypesPAI()

const std::vector< G4String > & G4EmParameters::TypesPAI ( ) const

Definition at line 1233 of file G4EmParameters.cc.

1234{
1235 return fBParameters->TypesPAI();
1236}

◆ TypesPhysics()

const std::vector< G4String > & G4EmParameters::TypesPhysics ( ) const

Definition at line 1277 of file G4EmParameters.cc.

1278{
1279 return fBParameters->TypesPhysics();
1280}

◆ Use3GammaAnnihilationOnFly()

G4bool G4EmParameters::Use3GammaAnnihilationOnFly ( ) const

◆ UseAngularGeneratorForIonisation()

G4bool G4EmParameters::UseAngularGeneratorForIonisation ( ) const

Definition at line 371 of file G4EmParameters.cc.

372{
373 return useAngGeneratorForIonisation;
374}

◆ UseCutAsFinalRange()

G4bool G4EmParameters::UseCutAsFinalRange ( ) const

Definition at line 231 of file G4EmParameters.cc.

232{
233 return cutAsFinalRange;
234}

◆ UseEPICS2017XS()

G4bool G4EmParameters::UseEPICS2017XS ( ) const

Definition at line 543 of file G4EmParameters.cc.

544{
545 return fUseEPICS2017XS;
546}

Referenced by G4BetheHeitlerModel::Initialise().

◆ UseICRU90Data()

G4bool G4EmParameters::UseICRU90Data ( ) const

Definition at line 427 of file G4EmParameters.cc.

428{
429 return fICRU90;
430}

Referenced by G4IonICRU73Data::Initialise(), and G4IonParametrisedLossModel::Initialise().

◆ UseMottCorrection()

G4bool G4EmParameters::UseMottCorrection ( ) const

Definition at line 382 of file G4EmParameters.cc.

383{
384 return useMottCorrection;
385}

Referenced by G4EmStandardPhysicsSS::ConstructProcess().

◆ UseRiGePairProductionModel()

G4bool G4EmParameters::UseRiGePairProductionModel ( ) const

Definition at line 565 of file G4EmParameters.cc.

566{
567 return fUseRiGePairProductionModel;
568}

Referenced by G4MuPairProduction::InitialiseEnergyLossProcess().

◆ Verbose()

◆ WorkerVerbose()

G4int G4EmParameters::WorkerVerbose ( ) const

Definition at line 1060 of file G4EmParameters.cc.

1061{
1062 return workerVerbose;
1063}

Referenced by G4DNABornIonisationModel::Initialise(), G4EmMultiModel::Initialise(), and G4EmTableUtil::PrepareEmProcess().

◆ operator<<

std::ostream & operator<< ( std::ostream & os,
const G4EmParameters & par )
friend

Definition at line 1603 of file G4EmParameters.cc.

1604{
1605 par.StreamInfo(os);
1606 return os;
1607}

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