Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmParameters Class Reference

#include <G4EmParameters.hh>

Public Member Functions

 ~G4EmParameters ()
 
void SetDefaults ()
 
std::ostream & StreamInfo (std::ostream &os) const
 
void Dump () const
 
void SetLossFluctuations (G4bool val)
 
G4bool LossFluctuation () const
 
void SetBuildCSDARange (G4bool val)
 
G4bool BuildCSDARange () const
 
void SetLPM (G4bool val)
 
G4bool LPM () const
 
void SetSpline (G4bool val)
 
G4bool Spline () const
 
void SetUseCutAsFinalRange (G4bool val)
 
G4bool UseCutAsFinalRange () const
 
void SetApplyCuts (G4bool val)
 
G4bool ApplyCuts () const
 
void SetFluo (G4bool val)
 
G4bool Fluo () const
 
void SetBeardenFluoDir (G4bool val)
 
G4bool BeardenFluoDir () const
 
void SetAuger (G4bool val)
 
G4bool Auger () const
 
void SetAugerCascade (G4bool val)
 
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 SetMuHadLateralDisplacement (G4bool val)
 
G4bool MuHadLateralDisplacement () const
 
void SetLatDisplacementBeyondSafety (G4bool val)
 
G4bool LatDisplacementBeyondSafety () 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 SetEmSaturation (G4EmSaturation *)
 
G4EmSaturationGetEmSaturation ()
 
void SetMinSubRange (G4double val)
 
G4double MinSubRange () 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 SetLinearLossLimit (G4double val)
 
G4double LinearLossLimit () const
 
void SetBremsstrahlungTh (G4double val)
 
G4double BremsstrahlungTh () const
 
void SetLambdaFactor (G4double val)
 
G4double LambdaFactor () const
 
void SetFactorForAngleLimit (G4double val)
 
G4double FactorForAngleLimit () const
 
void SetMscThetaLimit (G4double val)
 
G4double MscThetaLimit () const
 
void SetMscRangeFactor (G4double val)
 
G4double MscRangeFactor () const
 
void SetMscMuHadRangeFactor (G4double val)
 
G4double MscMuHadRangeFactor () const
 
void SetMscGeomFactor (G4double val)
 
G4double MscGeomFactor () const
 
void SetMscSkin (G4double val)
 
G4double MscSkin () const
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetStepFunctionMuHad (G4double v1, G4double v2)
 
void SetNumberOfBins (G4int val)
 
G4int NumberOfBins () const
 
void SetNumberOfBinsPerDecade (G4int val)
 
G4int NumberOfBinsPerDecade () const
 
void SetVerbose (G4int val)
 
G4int Verbose () const
 
void SetWorkerVerbose (G4int val)
 
G4int WorkerVerbose () const
 
void SetMscStepLimitType (G4MscStepLimitType val)
 
G4MscStepLimitType MscStepLimitType () const
 
void SetMscMuHadStepLimitType (G4MscStepLimitType val)
 
G4MscStepLimitType MscMuHadStepLimitType () const
 
void SetNuclearFormfactorType (G4NuclearFormfactorType val)
 
G4NuclearFormfactorType NuclearFormfactorType () const
 
void SetPIXECrossSectionModel (const G4String &)
 
const G4StringPIXECrossSectionModel ()
 
void SetPIXEElectronCrossSectionModel (const G4String &)
 
const G4StringPIXEElectronCrossSectionModel ()
 
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 AddMsc (const G4String &region, const G4String &type)
 
const std::vector< G4String > & RegionsMsc () const
 
const std::vector< G4String > & TypesMsc () const
 
void SetSubCutoff (G4bool val, 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 DefineRegParamForLoss (G4VEnergyLossProcess *, G4bool isElectron) const
 
void DefineRegParamForEM (G4VEmProcess *) const
 
void DefineRegParamForDeex (G4VAtomDeexcitation *) const
 

Static Public Member Functions

static G4EmParametersInstance ()
 

Friends

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

Detailed Description

Definition at line 72 of file G4EmParameters.hh.

Constructor & Destructor Documentation

G4EmParameters::~G4EmParameters ( )

Definition at line 78 of file G4EmParameters.cc.

79 {
80  delete theMessenger;
81  delete emSaturation;
82 }

Member Function Documentation

void G4EmParameters::ActivateAngularGeneratorForIonisation ( G4bool  val)

Definition at line 324 of file G4EmParameters.cc.

325 {
326  if(IsLocked()) { return; }
327  useAngGeneratorForIonisation = val;
328 }

Here is the caller graph for this function:

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

Definition at line 985 of file G4EmParameters.cc.

989 {
990  if(IsLocked()) { return; }
991  G4String r = CheckRegion(region);
992  if(length >= 0.0) {
993  G4int n = m_procForced.size();
994  for(G4int i=0; i<n; ++i) {
995  if(procname == m_procForced[i] && r == m_regnamesForced[i] ) {
996  m_lengthForced[i] = length;
997  m_weightForced[i]= wflag;
998  return;
999  }
1000  }
1001  m_regnamesForced.push_back(r);
1002  m_procForced.push_back(procname);
1003  m_lengthForced.push_back(length);
1004  m_weightForced.push_back(wflag);
1005  } else {
1007  ed << "Process: " << procname << " in region " << r
1008  << " : forced interacttion length= "
1009  << length << " is negative - ignored";
1010  PrintWarning(ed);
1011  }
1012 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

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

Definition at line 1015 of file G4EmParameters.cc.

1019 {
1020  if(IsLocked()) { return; }
1021  G4String r = CheckRegion(region);
1022  if(factor >= 0.0 && energyLimit >= 0.0) {
1023  G4int n = m_procBiasedSec.size();
1024  for(G4int i=0; i<n; ++i) {
1025  if(procname == m_procBiasedSec[i] && r == m_regnamesBiasedSec[i] ) {
1026  m_factBiasedSec[i] = factor;
1027  m_elimBiasedSec[i] = energyLimit;
1028  return;
1029  }
1030  }
1031  m_regnamesBiasedSec.push_back(r);
1032  m_procBiasedSec.push_back(procname);
1033  m_factBiasedSec.push_back(factor);
1034  m_elimBiasedSec.push_back(energyLimit);
1035  } else {
1037  ed << "Process: " << procname << " in region " << r
1038  << " : secondary bised factor= "
1039  << factor << ", Elim= " << energyLimit << " - ignored";
1040  PrintWarning(ed);
1041  }
1042 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

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

Definition at line 873 of file G4EmParameters.cc.

874 {
875  G4String r = CheckRegion(region);
876  G4int nreg = m_regnamesDNA.size();
877  for(G4int i=0; i<nreg; ++i) {
878  if(r == m_regnamesDNA[i]) { return; }
879  }
880  m_regnamesDNA.push_back(r);
881  m_typesDNA.push_back(type);
882 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

void G4EmParameters::AddMicroElec ( const G4String region)

Definition at line 858 of file G4EmParameters.cc.

859 {
860  G4String r = CheckRegion(region);
861  G4int nreg = m_regnamesME.size();
862  for(G4int i=0; i<nreg; ++i) {
863  if(r == m_regnamesME[i]) { return; }
864  }
865  m_regnamesME.push_back(r);
866 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

void G4EmParameters::AddMsc ( const G4String region,
const G4String type 
)

Definition at line 894 of file G4EmParameters.cc.

895 {
896  G4String r = CheckRegion(region);
897  G4int nreg = m_regnamesMsc.size();
898  for(G4int i=0; i<nreg; ++i) {
899  if(r == m_regnamesMsc[i]) { return; }
900  }
901  m_regnamesMsc.push_back(r);
902  m_typesMsc.push_back(type);
903 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

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

Definition at line 818 of file G4EmParameters.cc.

821 {
822  G4String r = CheckRegion(region);
823  G4int nreg = m_regnamesPAI.size();
824  for(G4int i=0; i<nreg; ++i) {
825  if((m_particlesPAI[i] == particle ||
826  m_particlesPAI[i] == "all" ||
827  particle == "all") &&
828  (m_regnamesPAI[i] == r ||
829  m_regnamesPAI[i] == "DefaultRegionForTheWorld" ||
830  r == "DefaultRegionForTheWorld") ) {
831 
832  m_typesPAI[i] = type;
833  if(particle == "all") { m_particlesPAI[i] = particle; }
834  if(r == "DefaultRegionForTheWorld") { m_regnamesPAI[i] = r; }
835  return;
836  }
837  }
838  m_particlesPAI.push_back(particle);
839  m_regnamesPAI.push_back(r);
840  m_typesPAI.push_back(type);
841 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

G4bool G4EmParameters::ApplyCuts ( ) const

Definition at line 217 of file G4EmParameters.cc.

218 {
219  return applyCuts;
220 }

Here is the caller graph for this function:

G4bool G4EmParameters::Auger ( ) const

Definition at line 251 of file G4EmParameters.cc.

252 {
253  return auger;
254 }

Here is the caller graph for this function:

G4bool G4EmParameters::AugerCascade ( ) const

Definition at line 263 of file G4EmParameters.cc.

264 {
265  return augerCascade;
266 }

Here is the caller graph for this function:

G4bool G4EmParameters::BeardenFluoDir ( ) const

Definition at line 239 of file G4EmParameters.cc.

240 {
241  return beardenFluoDir;
242 }

Here is the caller graph for this function:

G4bool G4EmParameters::BirksActive ( ) const

Definition at line 367 of file G4EmParameters.cc.

368 {
369  return birks;
370 }
G4double G4EmParameters::BremsstrahlungTh ( ) const

Definition at line 527 of file G4EmParameters.cc.

528 {
529  return bremsTh;
530 }

Here is the caller graph for this function:

G4bool G4EmParameters::BuildCSDARange ( ) const

Definition at line 173 of file G4EmParameters.cc.

174 {
175  return buildCSDARange;
176 }

Here is the caller graph for this function:

G4bool G4EmParameters::DeexcitationIgnoreCut ( ) const

Definition at line 286 of file G4EmParameters.cc.

287 {
288  return deexIgnoreCut;
289 }

Here is the caller graph for this function:

void G4EmParameters::DefineRegParamForDeex ( G4VAtomDeexcitation ptr) const

Definition at line 1114 of file G4EmParameters.cc.

1115 {
1116  G4int n = m_regnamesDeex.size();
1117  for(G4int i=0; i<n; ++i) {
1118  ptr->SetDeexcitationActiveRegion(m_regnamesDeex[i],
1119  m_fluo[i], m_auger[i], m_pixe[i]);
1120  }
1121 }
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::DefineRegParamForEM ( G4VEmProcess ptr) const

Definition at line 1084 of file G4EmParameters.cc.

1085 {
1086  G4int n = m_procBiasedXS.size();
1087  for(G4int i=0; i<n; ++i) {
1088  if(ptr->GetProcessName() == m_procBiasedXS[i]) {
1089  ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
1090  m_weightBiasedXS[i]);
1091  break;
1092  }
1093  }
1094  n = m_procForced.size();
1095  for(G4int i=0; i<n; ++i) {
1096  if(ptr->GetProcessName() == m_procForced[i]) {
1097  ptr->ActivateForcedInteraction(m_lengthForced[i],
1098  m_regnamesForced[i],
1099  m_weightForced[i]);
1100  break;
1101  }
1102  }
1103  n = m_procBiasedSec.size();
1104  for(G4int i=0; i<n; ++i) {
1105  if(ptr->GetProcessName() == m_procBiasedSec[i]) {
1106  ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
1107  m_factBiasedSec[i],
1108  m_elimBiasedSec[i]);
1109  break;
1110  }
1111  }
1112 }
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
int G4int
Definition: G4Types.hh:78
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::DefineRegParamForLoss ( G4VEnergyLossProcess ptr,
G4bool  isElectron 
) const

Definition at line 1044 of file G4EmParameters.cc.

1046 {
1047  if(isElectron) { ptr->SetStepFunction(dRoverRange, finalRange, false); }
1048  else { ptr->SetStepFunction(dRoverRangeMuHad, finalRangeMuHad, false); }
1049 
1050  G4RegionStore* regionStore = G4RegionStore::GetInstance();
1051  G4int n = m_regnamesSubCut.size();
1052  for(G4int i=0; i<n; ++i) {
1053  const G4Region* reg = regionStore->GetRegion(m_regnamesSubCut[i], false);
1054  if(reg) { ptr->ActivateSubCutoff(m_subCuts[i], reg); }
1055  }
1056  n = m_procBiasedXS.size();
1057  for(G4int i=0; i<n; ++i) {
1058  if(ptr->GetProcessName() == m_procBiasedXS[i]) {
1059  ptr->SetCrossSectionBiasingFactor(m_factBiasedXS[i],
1060  m_weightBiasedXS[i]);
1061  break;
1062  }
1063  }
1064  n = m_procForced.size();
1065  for(G4int i=0; i<n; ++i) {
1066  if(ptr->GetProcessName() == m_procForced[i]) {
1067  ptr->ActivateForcedInteraction(m_lengthForced[i],
1068  m_regnamesForced[i],
1069  m_weightForced[i]);
1070  break;
1071  }
1072  }
1073  n = m_procBiasedSec.size();
1074  for(G4int i=0; i<n; ++i) {
1075  if(ptr->GetProcessName() == m_procBiasedSec[i]) {
1076  ptr->ActivateSecondaryBiasing(m_regnamesBiasedSec[i],
1077  m_factBiasedSec[i],
1078  m_elimBiasedSec[i]);
1079  break;
1080  }
1081  }
1082 }
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4bool isElectron(G4int ityp)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
int G4int
Definition: G4Types.hh:78
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
static const G4double reg
static G4RegionStore * GetInstance()
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void ActivateSubCutoff(G4bool val, const G4Region *region=nullptr)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::Dump ( ) const

Definition at line 1199 of file G4EmParameters.cc.

1200 {
1201  StreamInfo(G4cout);
1202 }
std::ostream & StreamInfo(std::ostream &os) const
G4GLOB_DLL std::ostream G4cout

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmParameters::FactorForAngleLimit ( ) const

Definition at line 563 of file G4EmParameters.cc.

564 {
565  return factorForAngleLimit;
566 }

Here is the caller graph for this function:

G4bool G4EmParameters::Fluo ( ) const

Definition at line 228 of file G4EmParameters.cc.

229 {
230  return fluo;
231 }

Here is the caller graph for this function:

G4EmSaturation * G4EmParameters::GetEmSaturation ( )

Definition at line 380 of file G4EmParameters.cc.

381 {
382  if(!emSaturation) { SetBirksActive(true); }
383  return emSaturation;
384 }
void SetBirksActive(G4bool val)

Here is the call graph for this function:

Here is the caller graph for this function:

G4EmParameters * G4EmParameters::Instance ( void  )
static

Definition at line 67 of file G4EmParameters.cc.

68 {
69  if(nullptr == theInstance) {
70  static G4EmParameters manager;
71  theInstance = &manager;
72  }
73  return theInstance;
74 }
G4bool G4EmParameters::Integral ( ) const

Definition at line 352 of file G4EmParameters.cc.

353 {
354  return integral;
355 }

Here is the caller graph for this function:

G4double G4EmParameters::LambdaFactor ( ) const

Definition at line 545 of file G4EmParameters.cc.

546 {
547  return lambdaFactor;
548 }

Here is the caller graph for this function:

G4bool G4EmParameters::LatDisplacementBeyondSafety ( ) const

Definition at line 319 of file G4EmParameters.cc.

320 {
321  return latDisplacementBeyondSafety;
322 }

Here is the caller graph for this function:

G4bool G4EmParameters::LateralDisplacement ( ) const

Definition at line 297 of file G4EmParameters.cc.

298 {
299  return lateralDisplacement;
300 }

Here is the caller graph for this function:

G4double G4EmParameters::LinearLossLimit ( ) const

Definition at line 509 of file G4EmParameters.cc.

510 {
511  return linLossLimit;
512 }

Here is the caller graph for this function:

G4bool G4EmParameters::LossFluctuation ( ) const

Definition at line 162 of file G4EmParameters.cc.

163 {
164  return lossFluctuation;
165 }

Here is the caller graph for this function:

G4double G4EmParameters::LowestElectronEnergy ( ) const

Definition at line 473 of file G4EmParameters.cc.

474 {
475  return lowestElectronEnergy;
476 }

Here is the caller graph for this function:

G4double G4EmParameters::LowestMuHadEnergy ( ) const

Definition at line 491 of file G4EmParameters.cc.

492 {
493  return lowestMuHadEnergy;
494 }

Here is the caller graph for this function:

G4bool G4EmParameters::LPM ( ) const

Definition at line 184 of file G4EmParameters.cc.

185 {
186  return flagLPM;
187 }

Here is the caller graph for this function:

G4double G4EmParameters::MaxEnergyForCSDARange ( ) const

Definition at line 455 of file G4EmParameters.cc.

456 {
457  return maxKinEnergyCSDA;
458 }

Here is the caller graph for this function:

G4double G4EmParameters::MaxKinEnergy ( ) const

Definition at line 437 of file G4EmParameters.cc.

438 {
439  return maxKinEnergy;
440 }

Here is the caller graph for this function:

G4double G4EmParameters::MinKinEnergy ( ) const

Definition at line 418 of file G4EmParameters.cc.

419 {
420  return minKinEnergy;
421 }

Here is the caller graph for this function:

G4double G4EmParameters::MinSubRange ( ) const

Definition at line 399 of file G4EmParameters.cc.

400 {
401  return minSubRange;
402 }

Here is the caller graph for this function:

G4double G4EmParameters::MscGeomFactor ( ) const

Definition at line 635 of file G4EmParameters.cc.

636 {
637  return geomFactor;
638 }

Here is the caller graph for this function:

G4double G4EmParameters::MscMuHadRangeFactor ( ) const

Definition at line 617 of file G4EmParameters.cc.

618 {
619  return rangeFactorMuHad;
620 }

Here is the caller graph for this function:

G4MscStepLimitType G4EmParameters::MscMuHadStepLimitType ( ) const

Definition at line 764 of file G4EmParameters.cc.

765 {
766  return mscStepLimitMuHad;
767 }

Here is the caller graph for this function:

G4double G4EmParameters::MscRangeFactor ( ) const

Definition at line 599 of file G4EmParameters.cc.

600 {
601  return rangeFactor;
602 }

Here is the caller graph for this function:

G4double G4EmParameters::MscSkin ( ) const

Definition at line 653 of file G4EmParameters.cc.

654 {
655  return skin;
656 }

Here is the caller graph for this function:

G4MscStepLimitType G4EmParameters::MscStepLimitType ( ) const

Definition at line 753 of file G4EmParameters.cc.

754 {
755  return mscStepLimit;
756 }

Here is the caller graph for this function:

G4double G4EmParameters::MscThetaLimit ( ) const

Definition at line 581 of file G4EmParameters.cc.

582 {
583  return thetaLimit;
584 }

Here is the caller graph for this function:

G4bool G4EmParameters::MuHadLateralDisplacement ( ) const

Definition at line 308 of file G4EmParameters.cc.

309 {
310  return muhadLateralDisplacement;
311 }

Here is the caller graph for this function:

G4NuclearFormfactorType G4EmParameters::NuclearFormfactorType ( ) const

Definition at line 776 of file G4EmParameters.cc.

777 {
778  return nucFormfactor;
779 }

Here is the caller graph for this function:

G4int G4EmParameters::NumberOfBins ( ) const

Definition at line 700 of file G4EmParameters.cc.

701 {
702  return nbins;
703 }
G4int G4EmParameters::NumberOfBinsPerDecade ( ) const

Definition at line 719 of file G4EmParameters.cc.

720 {
721  return nbinsPerDecade;
722 }

Here is the caller graph for this function:

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

Definition at line 843 of file G4EmParameters.cc.

844 {
845  return m_particlesPAI;
846 }
G4bool G4EmParameters::Pixe ( ) const

Definition at line 275 of file G4EmParameters.cc.

276 {
277  return pixe;
278 }

Here is the caller graph for this function:

const G4String & G4EmParameters::PIXECrossSectionModel ( )

Definition at line 788 of file G4EmParameters.cc.

789 {
790  return namePIXE;
791 }

Here is the caller graph for this function:

const G4String & G4EmParameters::PIXEElectronCrossSectionModel ( )

Definition at line 799 of file G4EmParameters.cc.

800 {
801  return nameElectronPIXE;
802 }

Here is the caller graph for this function:

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

Definition at line 884 of file G4EmParameters.cc.

885 {
886  return m_regnamesDNA;
887 }

Here is the caller graph for this function:

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

Definition at line 868 of file G4EmParameters.cc.

869 {
870  return m_regnamesME;
871 }

Here is the caller graph for this function:

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

Definition at line 905 of file G4EmParameters.cc.

906 {
907  return m_regnamesMsc;
908 }

Here is the caller graph for this function:

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

Definition at line 848 of file G4EmParameters.cc.

849 {
850  return m_regnamesPAI;
851 }

Here is the caller graph for this function:

void G4EmParameters::SetApplyCuts ( G4bool  val)

Definition at line 211 of file G4EmParameters.cc.

212 {
213  if(IsLocked()) { return; }
214  applyCuts = val;
215 }

Here is the caller graph for this function:

void G4EmParameters::SetAuger ( G4bool  val)

Definition at line 244 of file G4EmParameters.cc.

245 {
246  if(IsLocked()) { return; }
247  auger = val;
248  if(val) { fluo = true; }
249 }

Here is the caller graph for this function:

void G4EmParameters::SetAugerCascade ( G4bool  val)

Definition at line 256 of file G4EmParameters.cc.

257 {
258  if(IsLocked()) { return; }
259  augerCascade = val;
260  if(val) { fluo = true; auger = true; }
261 }

Here is the caller graph for this function:

void G4EmParameters::SetBeardenFluoDir ( G4bool  val)

Definition at line 233 of file G4EmParameters.cc.

234 {
235  if(IsLocked()) { return; }
236  beardenFluoDir = val;
237 }

Here is the caller graph for this function:

void G4EmParameters::SetBirksActive ( G4bool  val)

Definition at line 357 of file G4EmParameters.cc.

358 {
359  if(IsLocked()) { return; }
360  birks = val;
361  if(birks) {
362  if(!emSaturation) { emSaturation = new G4EmSaturation(1); }
363  emSaturation->InitialiseG4Saturation();
364  }
365 }
void InitialiseG4Saturation()

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetBremsstrahlungTh ( G4double  val)

Definition at line 514 of file G4EmParameters.cc.

515 {
516  if(IsLocked()) { return; }
517  if(val > 0.0) {
518  bremsTh = val;
519  } else {
521  ed << "Value of bremsstrahlung threshold is out of range: "
522  << val/GeV << " GeV is ignored";
523  PrintWarning(ed);
524  }
525 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the caller graph for this function:

void G4EmParameters::SetBuildCSDARange ( G4bool  val)

Definition at line 167 of file G4EmParameters.cc.

168 {
169  if(IsLocked()) { return; }
170  buildCSDARange = val;
171 }

Here is the caller graph for this function:

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

Definition at line 931 of file G4EmParameters.cc.

933 {
934  if(IsLocked()) { return; }
935  if(fdeex) { fluo = true; }
936  G4String r = CheckRegion(region);
937  G4int nreg = m_regnamesDeex.size();
938  if(0 == nreg && r != "DefaultRegionForTheWorld") {
939  m_regnamesDeex.push_back("DefaultRegionForTheWorld");
940  m_fluo.push_back(false);
941  m_auger.push_back(false);
942  m_pixe.push_back(false);
943  nreg = 1;
944  }
945  for(G4int i=0; i<nreg; ++i) {
946  if(r == m_regnamesDeex[i]) {
947  m_fluo[i] = fdeex;
948  m_auger[i]= fauger;
949  m_pixe[i] = fpixe;
950  return;
951  }
952  }
953  m_regnamesDeex.push_back(r);
954  m_fluo.push_back(fdeex);
955  m_auger.push_back(fauger);
956  m_pixe.push_back(fpixe);
957 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

void G4EmParameters::SetDeexcitationIgnoreCut ( G4bool  val)

Definition at line 280 of file G4EmParameters.cc.

281 {
282  if(IsLocked()) { return; }
283  deexIgnoreCut = val;
284 }

Here is the caller graph for this function:

void G4EmParameters::SetDefaults ( )

Definition at line 96 of file G4EmParameters.cc.

97 {
98  if(!IsLocked()) { Initialise(); }
99 }

Here is the caller graph for this function:

void G4EmParameters::SetEmSaturation ( G4EmSaturation ptr)

Definition at line 372 of file G4EmParameters.cc.

373 {
374  if(emSaturation != ptr) {
375  delete emSaturation;
376  emSaturation = ptr;
377  }
378 }
void G4EmParameters::SetFactorForAngleLimit ( G4double  val)

Definition at line 550 of file G4EmParameters.cc.

551 {
552  if(IsLocked()) { return; }
553  if(val > 0.0) {
554  factorForAngleLimit = val;
555  } else {
557  ed << "Value of factor for enegry limit is out of range: "
558  << val << " is ignored";
559  PrintWarning(ed);
560  }
561 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetFluo ( G4bool  val)

Definition at line 222 of file G4EmParameters.cc.

223 {
224  if(IsLocked()) { return; }
225  fluo = val;
226 }

Here is the caller graph for this function:

void G4EmParameters::SetIntegral ( G4bool  val)

Definition at line 346 of file G4EmParameters.cc.

347 {
348  if(IsLocked()) { return; }
349  integral = val;
350 }

Here is the caller graph for this function:

void G4EmParameters::SetLambdaFactor ( G4double  val)

Definition at line 532 of file G4EmParameters.cc.

533 {
534  if(IsLocked()) { return; }
535  if(val > 0.0 && val < 1.0) {
536  lambdaFactor = val;
537  } else {
539  ed << "Value of lambda factor is out of range: " << val
540  << " is ignored";
541  PrintWarning(ed);
542  }
543 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetLatDisplacementBeyondSafety ( G4bool  val)

Definition at line 313 of file G4EmParameters.cc.

314 {
315  if(IsLocked()) { return; }
316  latDisplacementBeyondSafety = val;
317 }

Here is the caller graph for this function:

void G4EmParameters::SetLateralDisplacement ( G4bool  val)

Definition at line 291 of file G4EmParameters.cc.

292 {
293  if(IsLocked()) { return; }
294  lateralDisplacement = val;
295 }

Here is the caller graph for this function:

void G4EmParameters::SetLinearLossLimit ( G4double  val)

Definition at line 496 of file G4EmParameters.cc.

497 {
498  if(IsLocked()) { return; }
499  if(val > 0.0 && val < 0.5) {
500  linLossLimit = val;
501  } else {
503  ed << "Value of linLossLimit is out of range: " << val
504  << " is ignored";
505  PrintWarning(ed);
506  }
507 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetLossFluctuations ( G4bool  val)

Definition at line 156 of file G4EmParameters.cc.

157 {
158  if(IsLocked()) { return; }
159  lossFluctuation = val;
160 }

Here is the caller graph for this function:

void G4EmParameters::SetLowestElectronEnergy ( G4double  val)

Definition at line 460 of file G4EmParameters.cc.

461 {
462  if(IsLocked()) { return; }
463  if(val >= 0.0) {
464  lowestElectronEnergy = val;
465  } else {
467  ed << "Value of lowestElectronEnergy is out of range: "
468  << val/MeV << " MeV is ignored";
469  PrintWarning(ed);
470  }
471 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the caller graph for this function:

void G4EmParameters::SetLowestMuHadEnergy ( G4double  val)

Definition at line 478 of file G4EmParameters.cc.

479 {
480  if(IsLocked()) { return; }
481  if(val >= 0.0) {
482  lowestMuHadEnergy = val;
483  } else {
485  ed << "Value of lowestMuHadEnergy is out of range: "
486  << val/MeV << " MeV is ignored";
487  PrintWarning(ed);
488  }
489 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the caller graph for this function:

void G4EmParameters::SetLPM ( G4bool  val)

Definition at line 178 of file G4EmParameters.cc.

179 {
180  if(IsLocked()) { return; }
181  flagLPM = val;
182 }

Here is the caller graph for this function:

void G4EmParameters::SetMaxEnergy ( G4double  val)

Definition at line 423 of file G4EmParameters.cc.

424 {
425  if(IsLocked()) { return; }
426  if(val > minKinEnergy && val < 1.e+7*TeV) {
427  maxKinEnergy = val;
428  nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
429  } else {
431  ed << "Value of MaxKinEnergy is out of range: "
432  << val/GeV << " GeV is ignored";
433  PrintWarning(ed);
434  }
435 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double TeV
Definition: G4SIunits.hh:218
int G4lrint(double ad)
Definition: templates.hh:163
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetMaxEnergyForCSDARange ( G4double  val)

Definition at line 442 of file G4EmParameters.cc.

443 {
444  if(IsLocked()) { return; }
445  if(val > minKinEnergy && val <= 100*TeV) {
446  maxKinEnergyCSDA = val;
447  } else {
449  ed << "Value of MaxKinEnergyCSDA is out of range: "
450  << val/GeV << " GeV is ignored";
451  PrintWarning(ed);
452  }
453 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double TeV
Definition: G4SIunits.hh:218
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the caller graph for this function:

void G4EmParameters::SetMinEnergy ( G4double  val)

Definition at line 404 of file G4EmParameters.cc.

405 {
406  if(IsLocked()) { return; }
407  if(val > 1.e-3*eV && val < maxKinEnergy) {
408  minKinEnergy = val;
409  nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
410  } else {
412  ed << "Value of MinKinEnergy is out of range: " << val/MeV
413  << " MeV is ignored";
414  PrintWarning(ed);
415  }
416 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double eV
Definition: G4SIunits.hh:215
int G4lrint(double ad)
Definition: templates.hh:163
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetMinSubRange ( G4double  val)

Definition at line 386 of file G4EmParameters.cc.

387 {
388  if(IsLocked()) { return; }
389  if(val > 0.0 && val < 1.0) {
390  minSubRange = val;
391  } else {
393  ed << "Value of MinSubRange is out of range (0 - 1): " << val
394  << " is ignored";
395  PrintWarning(ed);
396  }
397 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetMscGeomFactor ( G4double  val)

Definition at line 622 of file G4EmParameters.cc.

623 {
624  if(IsLocked()) { return; }
625  if(val >= 1.0) {
626  geomFactor = val;
627  } else {
629  ed << "Value of geomFactor is out of range: "
630  << val << " is ignored";
631  PrintWarning(ed);
632  }
633 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetMscMuHadRangeFactor ( G4double  val)

Definition at line 604 of file G4EmParameters.cc.

605 {
606  if(IsLocked()) { return; }
607  if(val > 0.0 && val < 1.0) {
608  rangeFactorMuHad = val;
609  } else {
611  ed << "Value of rangeFactorMuHad is out of range: "
612  << val << " is ignored";
613  PrintWarning(ed);
614  }
615 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetMscMuHadStepLimitType ( G4MscStepLimitType  val)

Definition at line 758 of file G4EmParameters.cc.

759 {
760  if(IsLocked()) { return; }
761  mscStepLimitMuHad = val;
762 }

Here is the caller graph for this function:

void G4EmParameters::SetMscRangeFactor ( G4double  val)

Definition at line 586 of file G4EmParameters.cc.

587 {
588  if(IsLocked()) { return; }
589  if(val > 0.0 && val < 1.0) {
590  rangeFactor = val;
591  } else {
593  ed << "Value of rangeFactor is out of range: "
594  << val << " is ignored";
595  PrintWarning(ed);
596  }
597 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetMscSkin ( G4double  val)

Definition at line 640 of file G4EmParameters.cc.

641 {
642  if(IsLocked()) { return; }
643  if(val >= 1.0) {
644  skin = val;
645  } else {
647  ed << "Value of skin is out of range: "
648  << val << " is ignored";
649  PrintWarning(ed);
650  }
651 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Here is the caller graph for this function:

void G4EmParameters::SetMscStepLimitType ( G4MscStepLimitType  val)

Definition at line 747 of file G4EmParameters.cc.

748 {
749  if(IsLocked()) { return; }
750  mscStepLimit = val;
751 }

Here is the caller graph for this function:

void G4EmParameters::SetMscThetaLimit ( G4double  val)

Definition at line 568 of file G4EmParameters.cc.

569 {
570  if(IsLocked()) { return; }
571  if(val >= 0.0 && val <= pi) {
572  thetaLimit = val;
573  } else {
575  ed << "Value of polar angle limit is out of range: "
576  << val << " is ignored";
577  PrintWarning(ed);
578  }
579 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double pi
Definition: G4SIunits.hh:75

Here is the caller graph for this function:

void G4EmParameters::SetMuHadLateralDisplacement ( G4bool  val)

Definition at line 302 of file G4EmParameters.cc.

303 {
304  if(IsLocked()) { return; }
305  muhadLateralDisplacement = val;
306 }

Here is the caller graph for this function:

void G4EmParameters::SetNuclearFormfactorType ( G4NuclearFormfactorType  val)

Definition at line 770 of file G4EmParameters.cc.

771 {
772  if(IsLocked()) { return; }
773  nucFormfactor = val;
774 }

Here is the caller graph for this function:

void G4EmParameters::SetNumberOfBins ( G4int  val)

Definition at line 686 of file G4EmParameters.cc.

687 {
688  if(IsLocked()) { return; }
689  if(val >= 5 && val < 10000000) {
690  nbins = val;
691  nbinsPerDecade = G4lrint(nbins/std::log10(maxKinEnergy/minKinEnergy));
692  } else {
694  ed << "Value of number of bins is out of range: "
695  << val << " is ignored";
696  PrintWarning(ed);
697  }
698 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetNumberOfBinsPerDecade ( G4int  val)

Definition at line 705 of file G4EmParameters.cc.

706 {
707  if(IsLocked()) { return; }
708  if(val >= 5 && val < 1000000) {
709  nbinsPerDecade = val;
710  nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
711  } else {
713  ed << "Value of number of bins per decade is out of range: "
714  << val << " is ignored";
715  PrintWarning(ed);
716  }
717 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetPixe ( G4bool  val)

Definition at line 268 of file G4EmParameters.cc.

269 {
270  if(IsLocked()) { return; }
271  pixe = val;
272  if(val) { fluo = true; }
273 }

Here is the caller graph for this function:

void G4EmParameters::SetPIXECrossSectionModel ( const G4String sss)

Definition at line 781 of file G4EmParameters.cc.

782 {
783  G4cout << "G4EmParameters::SetPIXECrossSectionModel " << sss << G4endl;
784  if(IsLocked()) { return; }
785  namePIXE = sss;
786 }
G4GLOB_DLL std::ostream G4cout
static const char sss[MAX_N_PAR+2]
Definition: Evaluator.cc:64
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

void G4EmParameters::SetPIXEElectronCrossSectionModel ( const G4String sss)

Definition at line 793 of file G4EmParameters.cc.

794 {
795  if(IsLocked()) { return; }
796  nameElectronPIXE = sss;
797 }
static const char sss[MAX_N_PAR+2]
Definition: Evaluator.cc:64

Here is the caller graph for this function:

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

Definition at line 960 of file G4EmParameters.cc.

962 {
963  if(IsLocked()) { return; }
964  if(val > 0.0) {
965  G4int n = m_procBiasedXS.size();
966  for(G4int i=0; i<n; ++i) {
967  if(procname == m_procBiasedXS[i]) {
968  m_factBiasedXS[i] = val;
969  m_weightBiasedXS[i]= wflag;
970  return;
971  }
972  }
973  m_procBiasedXS.push_back(procname);
974  m_factBiasedXS.push_back(val);
975  m_weightBiasedXS.push_back(wflag);
976  } else {
978  ed << "Process: " << procname << " XS biasing factor "
979  << val << " is negative - ignored";
980  PrintWarning(ed);
981  }
982 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

void G4EmParameters::SetSpline ( G4bool  val)

Definition at line 189 of file G4EmParameters.cc.

190 {
191  if(IsLocked()) { return; }
192  spline = val;
193 }

Here is the caller graph for this function:

void G4EmParameters::SetStepFunction ( G4double  v1,
G4double  v2 
)

Definition at line 658 of file G4EmParameters.cc.

659 {
660  if(IsLocked()) { return; }
661  if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
662  dRoverRange = v1;
663  finalRange = v2;
664  } else {
666  ed << "Values of step function are out of range: "
667  << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
668  PrintWarning(ed);
669  }
670 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double mm
Definition: SystemOfUnits.h:95

Here is the caller graph for this function:

void G4EmParameters::SetStepFunctionMuHad ( G4double  v1,
G4double  v2 
)

Definition at line 672 of file G4EmParameters.cc.

673 {
674  if(IsLocked()) { return; }
675  if(v1 > 0.0 && v1 <= 1.0 && v2 > 0.0) {
676  dRoverRangeMuHad = v1;
677  finalRangeMuHad = v2;
678  } else {
680  ed << "Values of step function are out of range: "
681  << v1 << ", " << v2/CLHEP::mm << " mm - are ignored";
682  PrintWarning(ed);
683  }
684 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static constexpr double mm
Definition: SystemOfUnits.h:95

Here is the caller graph for this function:

void G4EmParameters::SetSubCutoff ( G4bool  val,
const G4String region = "" 
)

Definition at line 915 of file G4EmParameters.cc.

916 {
917  if(IsLocked()) { return; }
918  G4String r = CheckRegion(region);
919  G4int nreg = m_regnamesSubCut.size();
920  for(G4int i=0; i<nreg; ++i) {
921  if(r == m_regnamesSubCut[i]) {
922  m_subCuts[i] = val;
923  return;
924  }
925  }
926  m_regnamesSubCut.push_back(r);
927  m_subCuts.push_back(val);
928 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

void G4EmParameters::SetUseCutAsFinalRange ( G4bool  val)

Definition at line 200 of file G4EmParameters.cc.

201 {
202  if(IsLocked()) { return; }
203  cutAsFinalRange = val;
204 }

Here is the caller graph for this function:

void G4EmParameters::SetUseMottCorrection ( G4bool  val)

Definition at line 335 of file G4EmParameters.cc.

336 {
337  if(IsLocked()) { return; }
338  useMottCorrection = val;
339 }

Here is the caller graph for this function:

void G4EmParameters::SetVerbose ( G4int  val)

Definition at line 724 of file G4EmParameters.cc.

725 {
726  if(IsLocked()) { return; }
727  verbose = val;
728  workerVerbose = std::min(workerVerbose, verbose);
729 }
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmParameters::SetWorkerVerbose ( G4int  val)

Definition at line 736 of file G4EmParameters.cc.

737 {
738  if(IsLocked()) { return; }
739  workerVerbose = val;
740 }

Here is the caller graph for this function:

G4bool G4EmParameters::Spline ( ) const

Definition at line 195 of file G4EmParameters.cc.

196 {
197  return spline;
198 }

Here is the caller graph for this function:

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

Definition at line 1123 of file G4EmParameters.cc.

1124 {
1125  G4int prec = os.precision(5);
1126  os << "=======================================================================" << "\n";
1127  os << "====== Electromagnetic Physics Parameters ========" << "\n";
1128  os << "=======================================================================" << "\n";
1129  os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1130  os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1131  os << "LPM effect enabled " <<flagLPM << "\n";
1132  os << "Spline of EM tables enabled " <<spline << "\n";
1133  os << "Use cut as a final range enabled " <<finalRange << "\n";
1134  os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1135  os << "Fluorescence enabled " <<fluo << "\n";
1136  os << "Fluorescence Bearden data files enabled " <<beardenFluoDir << "\n";
1137  os << "Auger electron production enabled " <<auger << "\n";
1138  os << "Auger cascade enabled " <<augerCascade << "\n";
1139  os << "PIXE atomic de-excitation enabled " <<pixe << "\n";
1140  os << "De-excitation module ignores cuts " <<deexIgnoreCut << "\n";
1141  os << "Msc lateraral displacement for e+- enabled " <<lateralDisplacement << "\n";
1142  os << "Msc lateraral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1143  os << "Msc lateraral displacement beyond geometry safety " <<latDisplacementBeyondSafety << "\n";
1144  os << "Enable angular generator interface "
1145  <<useAngGeneratorForIonisation << "\n";
1146  os << "Use Mott correction for e- scattering "
1147  <<useMottCorrection << "\n";
1148  os << "Use integral approach for tracking "
1149  <<integral << "\n";
1150  os << "Use built-in Birks satuaration "
1151  << birks << "\n";
1152 
1153  os << "Factor of cut reduction for sub-cutoff method " <<minSubRange << "\n";
1154  os << "Min kinetic energy for tables "
1155  <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1156  os << "Max kinetic energy for tables "
1157  <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1158  os << "Max kinetic energy for CSDA tables "
1159  <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1160  os << "Lowest e+e- kinetic energy "
1161  <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1162  os << "Lowest muon/hadron kinetic energy "
1163  <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1164  os << "Linear loss limit " <<linLossLimit << "\n";
1165  os << "Bremsstrahlung energy threshold above which \n"
1166  << " primary is added to the list of secondary "
1167  <<G4BestUnit(bremsTh,"Energy") << "\n";
1168  os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1169  os << "Factor used for dynamic computation of angular \n"
1170  << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1171  os << "Fixed angular limit between single \n"
1172  << " and multiple scattering "
1173  <<thetaLimit/rad << " rad" << "\n";
1174  os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1175  os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1176  os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1177  os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1178  os << "Step function for e+- " <<"("<< dRoverRange
1179  << ", " << finalRange << " mm)\n";
1180  os << "Step function for muons/hadrons " <<"("<< dRoverRangeMuHad
1181  << ", " << finalRangeMuHad << " mm)\n";
1182 
1183  os << "Number of bins in tables " <<nbins << "\n";
1184  os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1185  os << "Verbose level " <<verbose << "\n";
1186  os << "Verbose level for worker thread " <<workerVerbose << "\n";
1187 
1188  os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1189  os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1190  os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1191 
1192  os << "Type of PIXE cross section for hadrons " <<namePIXE << "\n";
1193  os << "Type of PIXE cross section for e+- " <<nameElectronPIXE << "\n";
1194  os << "=======================================================================" << "\n";
1195  os.precision(prec);
1196  return os;
1197 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static constexpr double rad
Definition: G4SIunits.hh:149
int G4int
Definition: G4Types.hh:78
static const double prec
Definition: RanecuEngine.cc:58

Here is the caller graph for this function:

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

Definition at line 889 of file G4EmParameters.cc.

890 {
891  return m_typesDNA;
892 }

Here is the caller graph for this function:

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

Definition at line 910 of file G4EmParameters.cc.

911 {
912  return m_typesMsc;
913 }
const std::vector< G4String > & G4EmParameters::TypesPAI ( ) const

Definition at line 853 of file G4EmParameters.cc.

854 {
855  return m_typesPAI;
856 }
G4bool G4EmParameters::UseAngularGeneratorForIonisation ( ) const

Definition at line 330 of file G4EmParameters.cc.

331 {
332  return useAngGeneratorForIonisation;
333 }

Here is the caller graph for this function:

G4bool G4EmParameters::UseCutAsFinalRange ( ) const

Definition at line 206 of file G4EmParameters.cc.

207 {
208  return cutAsFinalRange;
209 }

Here is the caller graph for this function:

G4bool G4EmParameters::UseMottCorrection ( ) const

Definition at line 341 of file G4EmParameters.cc.

342 {
343  return useMottCorrection;
344 }
G4int G4EmParameters::Verbose ( ) const

Definition at line 731 of file G4EmParameters.cc.

732 {
733  return verbose;
734 }

Here is the caller graph for this function:

G4int G4EmParameters::WorkerVerbose ( ) const

Definition at line 742 of file G4EmParameters.cc.

743 {
744  return workerVerbose;
745 }

Here is the caller graph for this function:

Friends And Related Function Documentation

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

Definition at line 1204 of file G4EmParameters.cc.

1205 {
1206  return par.StreamInfo(os);
1207 }
std::ostream & StreamInfo(std::ostream &os) const

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