Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEmModel Class Referenceabstract

#include <G4VEmModel.hh>

+ Inheritance diagram for G4VEmModel:

Public Member Functions

 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 102 of file G4VEmModel.hh.

Constructor & Destructor Documentation

◆ G4VEmModel() [1/2]

G4VEmModel::G4VEmModel ( const G4String nam)
explicit

Definition at line 65 of file G4VEmModel.cc.

65 :
66 flucModel(nullptr),anglModel(nullptr), name(nam), lowLimit(0.1*CLHEP::keV),
67 highLimit(100.0*CLHEP::TeV),eMinActive(0.0),eMaxActive(DBL_MAX),
68 polarAngleLimit(CLHEP::pi),secondaryThreshold(DBL_MAX),
69 theLPMflag(false),flagDeexcitation(false),flagForceBuildTable(false),
70 isMaster(true),fElementData(nullptr),pParticleChange(nullptr),
71 xSectionTable(nullptr),pBaseMaterial(nullptr),idxTable(0),
72 lossFlucFlag(true),inveplus(1.0/CLHEP::eplus),pFactor(1.0),
73 fCurrentCouple(nullptr),fCurrentElement(nullptr),fCurrentIsotope(nullptr),
74 fTripletModel(nullptr),nsec(5)
75{
76 xsec.resize(nsec);
77 nSelectors = 0;
78 elmSelectors = nullptr;
79 localElmSelectors = true;
80 localTable = true;
81 useAngularGenerator = false;
82 useBaseMaterials = true;
83 isLocked = false;
84
85 fEmManager = G4LossTableManager::Instance();
86 fEmManager->Register(this);
87 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
90}
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void Register(G4VEnergyLossProcess *p)
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:440
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:442
size_t idxTable
Definition: G4VEmModel.hh:444
G4double inveplus
Definition: G4VEmModel.hh:446
G4bool lossFlucFlag
Definition: G4VEmModel.hh:445
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:439
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:441
G4double pFactor
Definition: G4VEmModel.hh:447
G4ElementData * fElementData
Definition: G4VEmModel.hh:438
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:443
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4VEmModel()

G4VEmModel::~G4VEmModel ( )
virtual

Definition at line 94 of file G4VEmModel.cc.

95{
96 if(localElmSelectors) {
97 for(G4int i=0; i<nSelectors; ++i) {
98 delete (*elmSelectors)[i];
99 }
100 delete elmSelectors;
101 }
102 delete anglModel;
103
104 if(localTable && xSectionTable != nullptr) {
106 delete xSectionTable;
107 xSectionTable = nullptr;
108 }
109 if(isMaster && fElementData != nullptr) {
110 delete fElementData;
111 fElementData = nullptr;
112 }
113 fEmManager->DeRegister(this);
114}
int G4int
Definition: G4Types.hh:85
void DeRegister(G4VEnergyLossProcess *p)
void clearAndDestroy()

◆ G4VEmModel() [2/2]

G4VEmModel::G4VEmModel ( const G4VEmModel )
delete

Member Function Documentation

◆ ChargeSquareRatio()

G4double G4VEmModel::ChargeSquareRatio ( const G4Track track)
virtual

Reimplemented in G4BraggIonGasModel, and G4BetheBlochIonGasModel.

Definition at line 383 of file G4VEmModel.cc.

384{
386 track.GetMaterial(), track.GetKineticEnergy());
387}
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391

Referenced by G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength().

◆ ComputeCrossSectionPerAtom() [1/2]

G4double G4VEmModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition part,
const G4Element elm,
G4double  kinEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inline

Definition at line 556 of file G4VEmModel.hh.

561{
563 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
564 cutEnergy,maxEnergy);
565}
G4double GetZ() const
Definition: G4Element.hh:130
G4double GetN() const
Definition: G4Element.hh:134
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:359
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:487

◆ ComputeCrossSectionPerAtom() [2/2]

G4double G4VEmModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented in G4eDPWACoulombScatteringModel, G4eBremsstrahlungRelModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4eSingleCoulombScatteringModel, G4IonCoulombScatteringModel, G4PolarizedComptonModel, G4eCoulombScatteringModel, G4hCoulombScatteringModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4PEEffectFluoModel, G4JAEAElasticScatteringModel, G4JAEAPolarizedElasticScatteringModel, G4LivermoreComptonModel, G4LivermoreComptonModifiedModel, G4LivermoreGammaConversionModelRC, G4LivermoreIonisationModel, G4LivermoreNuclearGammaConversionModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedGammaConversionModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermorePolarizedRayleighModel, G4LivermoreRayleighModel, G4LowEPComptonModel, G4LowEPPolarizedComptonModel, G4PenelopeAnnihilationModel, G4PenelopeGammaConversionModel, G4PenelopePhotoElectricModel, G4PenelopeRayleighModel, G4XrayRayleighModel, G4PenelopeRayleighModelMI, G4BoldyshevTripletModel, G4BetheHeitlerModel, G4PairProductionRelModel, G4eplusTo2GammaOKVIModel, G4eplusTo3GammaOKVIModel, G4eeToTwoGammaModel, G4EmMultiModel, G4LivermoreGammaConversion5DModel, G4LivermoreGammaConversionModel, G4WentzelVIModel, G4WentzelVIRelModel, G4mplIonisationWithDeltaModel, G4MuBetheBlochModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4AtimaEnergyLossModel, G4BetheBlochModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4LindhardSorensenIonModel, G4MollerBhabhaModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4eBremParametrizedModel, G4IonParametrisedLossModel, G4PenelopeComptonModel, G4PenelopeIonisationModel, G4UrbanAdjointMscModel, G4UrbanMscModel, and G4PenelopeBremsstrahlungModel.

Definition at line 359 of file G4VEmModel.cc.

362{
363 return 0.0;
364}

Referenced by ComputeCrossSectionPerAtom(), G4EmCalculator::ComputeCrossSectionPerAtom(), G4VEmProcess::ComputeCrossSectionPerAtom(), G4EmCalculator::ComputeCrossSectionPerShell(), G4AdjointIonIonisationModel::CorrectPostStepWeight(), CrossSectionPerVolume(), G4AdjointComptonModel::DiffCrossSectionPerAtomPrimToScatPrim(), G4AdjointhIonisationModel::DiffCrossSectionPerAtomPrimToSecond(), G4AdjointIonIonisationModel::DiffCrossSectionPerAtomPrimToSecond(), G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToSecond(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecondApproximated2(), and G4EmElementSelector::Initialise().

◆ ComputeCrossSectionPerShell()

G4double G4VEmModel::ComputeCrossSectionPerShell ( const G4ParticleDefinition ,
G4int  Z,
G4int  shellIdx,
G4double  kinEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Definition at line 369 of file G4VEmModel.cc.

372{
373 return 0.0;
374}

Referenced by G4EmCalculator::ComputeCrossSectionPerShell().

◆ ComputeDEDX()

G4double G4VEmModel::ComputeDEDX ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = DBL_MAX 
)
inlinevirtual

Reimplemented in G4EmMultiModel.

Definition at line 518 of file G4VEmModel.hh.

522{
523 SetCurrentCouple(couple);
524 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
525}
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245

◆ ComputeDEDXPerVolume()

◆ ComputeMeanFreePath()

G4double G4VEmModel::ComputeMeanFreePath ( const G4ParticleDefinition part,
G4double  kineticEnergy,
const G4Material material,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inline

Definition at line 543 of file G4VEmModel.hh.

548{
549 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
550 return (cross > 0.0) ? 1./cross : DBL_MAX;
551}
double G4double
Definition: G4Types.hh:83
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254

◆ CorrectionsAlongStep()

◆ CrossSection()

G4double G4VEmModel::CrossSection ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inline

Definition at line 529 of file G4VEmModel.hh.

534{
535 SetCurrentCouple(couple);
536 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
537 cutEnergy,maxEnergy);
538}

Referenced by G4EmModelManager::FillLambdaVector(), and G4EmMultiModel::SampleSecondaries().

◆ CrossSectionPerVolume()

G4double G4VEmModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented in G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4PAIModel, G4PAIPhotModel, G4BetheBlochNoDeltaModel, G4BraggNoDeltaModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4ICRU73NoDeltaModel, G4MuBetheBlochModel, G4AtimaEnergyLossModel, G4BetheBlochModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4LindhardSorensenIonModel, G4MollerBhabhaModel, G4PEEffectFluoModel, G4PenelopeRayleighModelMI, G4LEPTSAttachmentModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSExcitationModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, G4LEPTSVibExcitationModel, G4PenelopeComptonModel, G4GoudsmitSaundersonMscModel, G4eplusTo2GammaOKVIModel, G4eplusTo3GammaOKVIModel, G4eeToTwoGammaModel, G4IonParametrisedLossModel, G4DNABornExcitationModel1, G4DNABornExcitationModel2, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNAChampionElasticModel, G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNADingfelderChargeDecreaseModel, G4DNADingfelderChargeIncreaseModel, G4DNAELSEPAElasticModel, G4DNAEmfietzoglouExcitationModel, G4DNAEmfietzoglouIonisationModel, G4DNAIonElasticModel, G4DNAMeltonAttachmentModel, G4DNAMillerGreenExcitationModel, G4DNAModelInterface, G4TDNAOneStepThermalizationModel< MODEL >, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNASancheExcitationModel, G4DNAScreenedRutherfordElasticModel, G4DNATransformElectronModel, G4DNAUeharaScreenedRutherfordElasticModel, G4MicroElecElasticModel, G4MicroElecInelasticModel, G4MicroElecElasticModel_new, G4MicroElecInelasticModel_new, G4MicroElecLOPhononModel, G4PenelopeBremsstrahlungModel, and G4PenelopeIonisationModel.

Definition at line 254 of file G4VEmModel.cc.

259{
260 SetupForMaterial(p, material, ekin);
261 G4double cross = 0.0;
262 const G4double* theAtomNumDensityVector =
263 material->GetVecNbOfAtomsPerVolume();
264 G4int nelm = material->GetNumberOfElements();
265 if(nelm > nsec) {
266 xsec.resize(nelm);
267 nsec = nelm;
268 }
269 for (G4int i=0; i<nelm; ++i) {
270 cross += theAtomNumDensityVector[i]*
271 ComputeCrossSectionPerAtom(p,material->GetElement(i),ekin,emin,emax);
272 xsec[i] = cross;
273 }
274 return cross;
275}
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:449

Referenced by G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4EmCalculator::ComputeCrossSectionPerVolume(), ComputeMeanFreePath(), CrossSection(), G4LivermorePhotoElectricModel::CrossSectionPerVolume(), G4LivermorePolarizedPhotoElectricModel::CrossSectionPerVolume(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4DNADummyModel::CrossSectionPerVolume(), G4VEmProcess::CrossSectionPerVolume(), G4VEnergyLossProcess::CrossSectionPerVolume(), G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToSecond(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecondApproximated1(), G4AdjointBremsstrahlungModel::GetAdjointCrossSection(), G4VMscModel::GetTransportMeanFreePath(), SelectRandomAtom(), and Value().

◆ CurrentCouple()

◆ DeexcitationFlag()

G4bool G4VEmModel::DeexcitationFlag ( ) const
inline

Definition at line 694 of file G4VEmModel.hh.

695{
696 return flagDeexcitation;
697}

Referenced by G4EmModelManager::DumpModelList().

◆ DefineForRegion()

void G4VEmModel::DefineForRegion ( const G4Region )
virtual

Reimplemented in G4PAIModel, and G4PAIPhotModel.

Definition at line 378 of file G4VEmModel.cc.

379{}

Referenced by G4EmModelManager::AddEmModel().

◆ ForceBuildTableFlag()

G4bool G4VEmModel::ForceBuildTableFlag ( ) const
inline

Definition at line 701 of file G4VEmModel.hh.

702{
703 return flagForceBuildTable;
704}

Referenced by G4VMscModel::GetParticleChangeForMSC().

◆ GetAngularDistribution()

G4VEmAngularDistribution * G4VEmModel::GetAngularDistribution ( )
inline

Definition at line 611 of file G4VEmModel.hh.

612{
613 return anglModel;
614}

Referenced by G4EmModelManager::DumpModelList(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), G4MollerBhabhaModel::Initialise(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4PAIModel::SampleSecondaries(), G4PAIPhotModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4AtimaEnergyLossModel::SampleSecondaries(), G4BetheBlochModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4BraggIonModel::SampleSecondaries(), G4BraggModel::SampleSecondaries(), G4ICRU73QOModel::SampleSecondaries(), G4LindhardSorensenIonModel::SampleSecondaries(), G4MollerBhabhaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), and G4IonParametrisedLossModel::SampleSecondaries().

◆ GetChargeSquareRatio()

◆ GetCrossSectionTable()

G4PhysicsTable * G4VEmModel::GetCrossSectionTable ( )
inline

◆ GetCurrentElement()

◆ GetCurrentIsotope()

const G4Isotope * G4VEmModel::GetCurrentIsotope ( ) const
inline

Definition at line 502 of file G4VEmModel.hh.

503{
504 return fCurrentIsotope;
505}

Referenced by G4VEmProcess::GetTargetIsotope().

◆ GetElementData()

G4ElementData * G4VEmModel::GetElementData ( )
inline

◆ GetElementSelectors()

◆ GetModelOfFluctuations()

◆ GetName()

◆ GetPartialCrossSection()

G4double G4VEmModel::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition ,
G4double  kineticEnergy 
)
virtual

◆ GetParticleChangeForGamma()

G4ParticleChangeForGamma * G4VEmModel::GetParticleChangeForGamma ( )
protected

Definition at line 133 of file G4VEmModel.cc.

134{
135 G4ParticleChangeForGamma* p = nullptr;
136 if (pParticleChange != nullptr) {
137 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
138 } else {
139 p = new G4ParticleChangeForGamma();
140 pParticleChange = p;
141 }
142 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
143 return p;
144}
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456

Referenced by G4MicroElecLOPhononModel::G4MicroElecLOPhononModel(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNATransformElectronModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4LEPTSAttachmentModel::Initialise(), G4LEPTSDissociationModel::Initialise(), G4LEPTSElasticModel::Initialise(), G4LEPTSExcitationModel::Initialise(), G4LEPTSIonisationModel::Initialise(), G4LEPTSPositroniumModel::Initialise(), G4LEPTSRotExcitationModel::Initialise(), G4LEPTSVibExcitationModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreComptonModifiedModel::Initialise(), G4LivermoreGammaConversionModelRC::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedPhotoElectricGDModel::Initialise(), G4LivermorePolarizedPhotoElectricModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4XrayRayleighModel::Initialise(), G4PolarizedAnnihilationModel::Initialise(), G4eplusTo3GammaOKVIModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4IonCoulombScatteringModel::Initialise(), G4eeToHadronsMultiModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4MicroElecLOPhononModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4BetheHeitlerModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4eeToTwoGammaModel::Initialise(), G4eplusTo2GammaOKVIModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4PairProductionRelModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAModelInterface::Initialise(), and G4DNAIonElasticModel::Initialise().

◆ GetParticleChangeForLoss()

◆ GetParticleCharge()

G4double G4VEmModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material ,
G4double  kineticEnergy 
)
virtual

◆ GetTripletModel()

G4VEmModel * G4VEmModel::GetTripletModel ( )
inline

Definition at line 628 of file G4VEmModel.hh.

629{
630 return fTripletModel;
631}

Referenced by G4eBremsstrahlungRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), and G4eBremsstrahlungRelModel::SampleSecondaries().

◆ HighEnergyActivationLimit()

G4double G4VEmModel::HighEnergyActivationLimit ( ) const
inline

◆ HighEnergyLimit()

G4double G4VEmModel::HighEnergyLimit ( ) const
inline

Definition at line 645 of file G4VEmModel.hh.

646{
647 return highLimit;
648}

Referenced by G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNATransformElectronModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eeToHadronsModel::G4eeToHadronsModel(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), G4VMscModel::GetParticleChangeForMSC(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreComptonModifiedModel::Initialise(), G4LivermoreGammaConversionModelRC::Initialise(), G4LivermoreIonisationModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4WentzelVIModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4EmModelManager::Initialise(), G4DNADummyModel::Initialise(), G4DNAIonElasticModel::Initialise(), G4mplIonisation::InitialiseEnergyLossProcess(), G4alphaIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4MuBremsstrahlungModel::InitialiseLocal(), G4MuPairProductionModel::InitialiseLocal(), G4eBremsstrahlungRelModel::InitialiseLocal(), G4PairProductionRelModel::InitialiseLocal(), G4CoulombScattering::InitialiseProcess(), G4VEmProcess::PostStepDoIt(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4EmConfigurator::SetExtraEmModel(), G4mplIonisationModel::SetParticle(), G4mplIonisationWithDeltaModel::SetParticle(), and G4eBremsstrahlung::StreamProcessInfo().

◆ Initialise()

virtual void G4VEmModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
pure virtual

Implemented in G4DNAChampionElasticModel, G4DNACPA100ElasticModel, G4DNADingfelderChargeDecreaseModel, G4DNADingfelderChargeIncreaseModel, G4DNAELSEPAElasticModel, G4DNAMeltonAttachmentModel, G4DNAMillerGreenExcitationModel, G4TDNAOneStepThermalizationModel< MODEL >, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNASancheExcitationModel, G4DNAScreenedRutherfordElasticModel, G4DNATransformElectronModel, G4DNAUeharaScreenedRutherfordElasticModel, G4LEPTSAttachmentModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSExcitationModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, G4LEPTSVibExcitationModel, G4BoldyshevTripletModel, G4IonParametrisedLossModel, G4JAEAElasticScatteringModel, G4JAEAPolarizedElasticScatteringModel, G4LivermoreBremsstrahlungModel, G4LivermoreComptonModel, G4LivermoreComptonModifiedModel, G4LivermoreGammaConversionModelRC, G4LivermoreIonisationModel, G4LivermoreNuclearGammaConversionModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedGammaConversionModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermorePolarizedPhotoElectricModel, G4LivermorePolarizedRayleighModel, G4LivermoreRayleighModel, G4LowEPComptonModel, G4LowEPPolarizedComptonModel, G4MicroElecElasticModel, G4MicroElecInelasticModel, G4PenelopeAnnihilationModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4PenelopeGammaConversionModel, G4PenelopeIonisationModel, G4PenelopePhotoElectricModel, G4PenelopeRayleighModel, G4GoudsmitSaundersonMscModel, G4XrayRayleighModel, G4PolarizedAnnihilationModel, G4eplusTo3GammaOKVIModel, G4eSingleCoulombScatteringModel, G4ICRU49NuclearStoppingModel, G4IonCoulombScatteringModel, G4PAIModel, G4PAIPhotModel, G4DummyModel, G4EmMultiModel, G4UrbanAdjointMscModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4mplIonisationModel, G4mplIonisationWithDeltaModel, G4LivermoreGammaConversion5DModel, G4LivermoreGammaConversionModel, G4MicroElecElasticModel_new, G4MicroElecInelasticModel_new, G4MicroElecLOPhononModel, G4PenelopeRayleighModelMI, G4MuBetheBlochModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4ePolarizedBremsstrahlungModel, G4PolarizedGammaConversionModel, G4AtimaEnergyLossModel, G4BetheBlochModel, G4BetheHeitler5DModel, G4BetheHeitlerModel, G4BraggIonModel, G4BraggModel, G4eBremParametrizedModel, G4eBremsstrahlungRelModel, G4eCoulombScatteringModel, G4eDPWACoulombScatteringModel, G4eeToTwoGammaModel, G4eplusTo2GammaOKVIModel, G4hCoulombScatteringModel, G4ICRU73QOModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4LindhardSorensenIonModel, G4MollerBhabhaModel, G4PairProductionRelModel, G4PEEffectFluoModel, G4SeltzerBergerModel, G4UrbanMscModel, G4WentzelVIModel, G4WentzelVIRelModel, G4DNABornExcitationModel1, G4DNABornExcitationModel2, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNAEmfietzoglouExcitationModel, G4DNAEmfietzoglouIonisationModel, G4DNAModelInterface, G4DNAIonElasticModel, and G4PolarizedPEEffectModel.

Referenced by G4eBremsstrahlungRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4EmModelManager::Initialise(), G4DNADummyModel::Initialise(), and InitialiseElementSelectors().

◆ InitialiseElementSelectors()

void G4VEmModel::InitialiseElementSelectors ( const G4ParticleDefinition part,
const G4DataVector cuts 
)

Definition at line 148 of file G4VEmModel.cc.

150{
151 // using spline for element selectors should be investigated in details
152 // because small number of points may provide biased results
153 // large number of points requires significant increase of memory
154 G4bool spline = false;
155
156 //G4cout << "IES: for " << GetName() << " Emin(MeV)= " << lowLimit/MeV
157 // << " Emax(MeV)= " << highLimit/MeV << G4endl;
158
159 // two times less bins because probability functon is normalized
160 // so correspondingly is more smooth
161 if(highLimit <= lowLimit) { return; }
162
164
165 G4ProductionCutsTable* theCoupleTable=
167 G4int numOfCouples = theCoupleTable->GetTableSize();
168
169 // prepare vector
170 if(!elmSelectors) {
171 elmSelectors = new std::vector<G4EmElementSelector*>;
172 }
173 if(numOfCouples > nSelectors) {
174 for(G4int i=nSelectors; i<numOfCouples; ++i) {
175 elmSelectors->push_back(nullptr);
176 }
177 nSelectors = numOfCouples;
178 }
179
180 // initialise vector
181 for(G4int i=0; i<numOfCouples; ++i) {
182
183 // no need in element selectors for infionite cuts
184 if(cuts[i] == DBL_MAX) { continue; }
185
186 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
187 auto material = couple->GetMaterial();
188 SetCurrentCouple(couple);
189
190 // selector already exist check if should be deleted
191 G4bool create = true;
192 if((*elmSelectors)[i]) {
193 if(material == ((*elmSelectors)[i])->GetMaterial()) { create = false; }
194 else { delete (*elmSelectors)[i]; }
195 }
196 if(create) {
197 G4double emin = std::max(lowLimit,
198 MinPrimaryEnergy(material, part, cuts[i]));
199 G4double emax = std::max(highLimit, 10*emin);
200 static const G4double invlog106 = 1.0/(6*G4Log(10.));
201 G4int nbins = (G4int)(nbinsPerDec*G4Log(emax/emin)*invlog106);
202 nbins = std::max(nbins, 3);
203
204 (*elmSelectors)[i] = new G4EmElementSelector(this,material,nbins,
205 emin,emax,spline);
206 }
207 ((*elmSelectors)[i])->Initialise(part, cuts[i]);
208 /*
209 G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
210 << " " << part->GetParticleName()
211 << " for " << GetName() << " cut= " << cuts[i]
212 << " " << (*elmSelectors)[i] << G4endl;
213 ((*elmSelectors)[i])->Dump(part);
214 */
215 }
216}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
bool G4bool
Definition: G4Types.hh:86
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:424
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0

Referenced by G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreComptonModifiedModel::Initialise(), G4LivermoreGammaConversionModelRC::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4XrayRayleighModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4PairProductionRelModel::Initialise(), and G4SeltzerBergerModel::Initialise().

◆ InitialiseForElement()

◆ InitialiseForMaterial()

void G4VEmModel::InitialiseForMaterial ( const G4ParticleDefinition part,
const G4Material material 
)
virtual

Definition at line 226 of file G4VEmModel.cc.

228{
229 if(material != nullptr) {
230 size_t n = material->GetNumberOfElements();
231 for(size_t i=0; i<n; ++i) {
232 G4int Z = material->GetElement(i)->GetZasInt();
233 InitialiseForElement(part, Z);
234 }
235 }
236}
G4int GetZasInt() const
Definition: G4Element.hh:131
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:240

◆ InitialiseLocal()

◆ IsActive()

◆ IsLocked()

◆ IsMaster()

G4bool G4VEmModel::IsMaster ( ) const
inline

Definition at line 736 of file G4VEmModel.hh.

737{
738 return isMaster;
739}

Referenced by G4PenelopePhotoElectricModel::GetNumberOfShellXS(), G4VMscModel::GetParticleChangeForMSC(), G4BoldyshevTripletModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreBremsstrahlungModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversionModelRC::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedPhotoElectricGDModel::Initialise(), G4LivermorePolarizedPhotoElectricModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4mplIonisationModel::Initialise(), G4mplIonisationWithDeltaModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4BetheBlochModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4eeToTwoGammaModel::Initialise(), G4eplusTo2GammaOKVIModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4UrbanMscModel::Initialise(), G4WentzelVIModel::Initialise(), G4PenelopeBremsstrahlungModel::InitialiseLocal(), G4BetheHeitlerModel::~G4BetheHeitlerModel(), G4BoldyshevTripletModel::~G4BoldyshevTripletModel(), G4BraggIonModel::~G4BraggIonModel(), G4BraggModel::~G4BraggModel(), G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel(), G4eDPWACoulombScatteringModel::~G4eDPWACoulombScatteringModel(), G4GoudsmitSaundersonMscModel::~G4GoudsmitSaundersonMscModel(), G4JAEAElasticScatteringModel::~G4JAEAElasticScatteringModel(), G4JAEAPolarizedElasticScatteringModel::~G4JAEAPolarizedElasticScatteringModel(), G4LivermoreBremsstrahlungModel::~G4LivermoreBremsstrahlungModel(), G4LivermoreComptonModel::~G4LivermoreComptonModel(), G4LivermoreGammaConversion5DModel::~G4LivermoreGammaConversion5DModel(), G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel(), G4LivermoreGammaConversionModelRC::~G4LivermoreGammaConversionModelRC(), G4LivermoreNuclearGammaConversionModel::~G4LivermoreNuclearGammaConversionModel(), G4LivermorePhotoElectricModel::~G4LivermorePhotoElectricModel(), G4LivermorePolarizedComptonModel::~G4LivermorePolarizedComptonModel(), G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel(), G4LivermorePolarizedPhotoElectricGDModel::~G4LivermorePolarizedPhotoElectricGDModel(), G4LivermorePolarizedPhotoElectricModel::~G4LivermorePolarizedPhotoElectricModel(), G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel(), G4LivermoreRayleighModel::~G4LivermoreRayleighModel(), G4LowEPComptonModel::~G4LowEPComptonModel(), G4LowEPPolarizedComptonModel::~G4LowEPPolarizedComptonModel(), G4mplIonisationModel::~G4mplIonisationModel(), G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel(), G4PAIModel::~G4PAIModel(), G4PAIPhotModel::~G4PAIPhotModel(), G4PairProductionRelModel::~G4PairProductionRelModel(), G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel(), G4PenelopeGammaConversionModel::~G4PenelopeGammaConversionModel(), G4PenelopeIonisationModel::~G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel(), G4PenelopeRayleighModel::~G4PenelopeRayleighModel(), G4PenelopeRayleighModelMI::~G4PenelopeRayleighModelMI(), G4SeltzerBergerModel::~G4SeltzerBergerModel(), G4UrbanMscModel::~G4UrbanMscModel(), and G4WentzelVIModel::~G4WentzelVIModel().

◆ LowEnergyActivationLimit()

G4double G4VEmModel::LowEnergyActivationLimit ( ) const
inline

◆ LowEnergyLimit()

G4double G4VEmModel::LowEnergyLimit ( ) const
inline

Definition at line 652 of file G4VEmModel.hh.

653{
654 return lowLimit;
655}

Referenced by G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom(), G4KleinNishinaCompton::ComputeCrossSectionPerAtom(), G4KleinNishinaModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModifiedModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron(), G4EmCalculator::ComputeDEDX(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(), G4IonParametrisedLossModel::ComputeDEDXPerVolume(), G4IonParametrisedLossModel::CorrectionsAlongStep(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4PenelopeComptonModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eeToHadronsModel::G4eeToHadronsModel(), G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(), G4VMscModel::GetParticleChangeForMSC(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreComptonModifiedModel::Initialise(), G4LivermoreGammaConversionModelRC::Initialise(), G4LivermoreIonisationModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuPairProductionModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4WentzelVIModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4EmModelManager::Initialise(), G4DNADummyModel::Initialise(), G4DNAIonElasticModel::Initialise(), G4mplIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlungRelModel::InitialiseLocal(), G4PairProductionRelModel::InitialiseLocal(), G4CoulombScattering::InitialiseProcess(), G4VEmProcess::PostStepDoIt(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4EmConfigurator::SetExtraEmModel(), G4mplIonisationModel::SetParticle(), and G4mplIonisationWithDeltaModel::SetParticle().

◆ LPMFlag()

◆ MaxSecondaryEnergy()

G4double G4VEmModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kineticEnergy 
)
protectedvirtual

◆ MaxSecondaryKinEnergy()

G4double G4VEmModel::MaxSecondaryKinEnergy ( const G4DynamicParticle dynParticle)
inline

◆ MinEnergyCut()

◆ MinPrimaryEnergy()

◆ ModelDescription()

void G4VEmModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented in G4eeToHadronsMultiModel.

Definition at line 478 of file G4VEmModel.cc.

479{
480 outFile << "The description for this model has not been written yet.\n";
481}

◆ operator=()

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

◆ PolarAngleLimit()

G4double G4VEmModel::PolarAngleLimit ( ) const
inline

◆ SampleSecondaries()

virtual void G4VEmModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin = 0.0,
G4double  tmax = DBL_MAX 
)
pure virtual

Implemented in G4LivermoreBremsstrahlungModel, G4eBremParametrizedModel, G4eBremsstrahlungRelModel, G4SeltzerBergerModel, G4DNABornExcitationModel1, G4DNABornExcitationModel2, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNAChampionElasticModel, G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNADingfelderChargeDecreaseModel, G4DNADingfelderChargeIncreaseModel, G4DNAELSEPAElasticModel, G4DNAEmfietzoglouExcitationModel, G4DNAEmfietzoglouIonisationModel, G4DNAIonElasticModel, G4DNAMeltonAttachmentModel, G4DNAMillerGreenExcitationModel, G4TDNAOneStepThermalizationModel< MODEL >, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNASancheExcitationModel, G4DNAScreenedRutherfordElasticModel, G4DNATransformElectronModel, G4DNAUeharaScreenedRutherfordElasticModel, G4BoldyshevTripletModel, G4JAEAElasticScatteringModel, G4JAEAPolarizedElasticScatteringModel, G4LivermoreComptonModel, G4LivermoreComptonModifiedModel, G4LivermoreGammaConversionModelRC, G4LivermoreIonisationModel, G4LivermoreNuclearGammaConversionModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedGammaConversionModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermorePolarizedPhotoElectricModel, G4LivermorePolarizedRayleighModel, G4LivermoreRayleighModel, G4LowEPComptonModel, G4LowEPPolarizedComptonModel, G4MicroElecElasticModel, G4MicroElecInelasticModel, G4PenelopeAnnihilationModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4PenelopeGammaConversionModel, G4PenelopeIonisationModel, G4PenelopePhotoElectricModel, G4PenelopeRayleighModel, G4XrayRayleighModel, G4PolarizedAnnihilationModel, G4eSingleCoulombScatteringModel, G4HeatedKleinNishinaCompton, G4IonCoulombScatteringModel, G4PAIModel, G4PAIPhotModel, G4mplIonisationModel, G4mplIonisationWithDeltaModel, G4MicroElecElasticModel_new, G4MicroElecInelasticModel_new, G4MicroElecLOPhononModel, G4PenelopeRayleighModelMI, G4MuBetheBlochModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4ePolarizedBremsstrahlungModel, G4PolarizedComptonModel, G4PolarizedGammaConversionModel, G4PolarizedMollerBhabhaModel, G4PolarizedPEEffectModel, G4AtimaEnergyLossModel, G4BetheBlochModel, G4BetheHeitlerModel, G4BraggIonModel, G4BraggModel, G4eCoulombScatteringModel, G4eDPWACoulombScatteringModel, G4eeToTwoGammaModel, G4hCoulombScatteringModel, G4ICRU73QOModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4LindhardSorensenIonModel, G4MollerBhabhaModel, G4PairProductionRelModel, G4PEEffectFluoModel, G4DummyModel, G4EmMultiModel, G4VMscModel, G4eplusTo2GammaOKVIModel, G4eplusTo3GammaOKVIModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4LEPTSAttachmentModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSExcitationModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, G4LEPTSVibExcitationModel, G4IonParametrisedLossModel, G4ICRU49NuclearStoppingModel, G4DNAModelInterface, and G4BetheHeitler5DModel.

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4DNADummyModel::SampleSecondaries(), and G4VEnergyLossProcess::SampleSubCutSecondaries().

◆ SecondaryThreshold()

◆ SelectIsotopeNumber()

G4int G4VEmModel::SelectIsotopeNumber ( const G4Element elm)

Definition at line 337 of file G4VEmModel.cc.

338{
340 size_t ni = elm->GetNumberOfIsotopes();
341 fCurrentIsotope = elm->GetIsotope(0);
342 size_t idx = 0;
343 if(ni > 1) {
344 const G4double* ab = elm->GetRelativeAbundanceVector();
346 for(; idx<ni; ++idx) {
347 x -= ab[idx];
348 if (x <= 0.0) {
349 fCurrentIsotope = elm->GetIsotope(idx);
350 break;
351 }
352 }
353 }
354 return fCurrentIsotope->GetN();
355}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4int GetN() const
Definition: G4Isotope.hh:93

Referenced by G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ SelectRandomAtom() [1/2]

const G4Element * G4VEmModel::SelectRandomAtom ( const G4Material material,
const G4ParticleDefinition pd,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)

Definition at line 293 of file G4VEmModel.cc.

298{
299 size_t n = material->GetNumberOfElements();
300 fCurrentElement = material->GetElement(0);
301 if (n > 1) {
303 G4VEmModel::CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax);
304 for(size_t i=0; i<n; ++i) {
305 if (x <= xsec[i]) {
306 fCurrentElement = material->GetElement(i);
307 break;
308 }
309 }
310 }
311 return fCurrentElement;
312}

◆ SelectRandomAtom() [2/2]

const G4Element * G4VEmModel::SelectRandomAtom ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inline

Definition at line 570 of file G4VEmModel.hh.

575{
576 SetCurrentCouple(couple);
577 fCurrentElement = (nSelectors > 0) ?
578 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
579 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
580 fCurrentIsotope = nullptr;
581 return fCurrentElement;
582}
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570

Referenced by G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricGDModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), SelectRandomAtom(), and SelectTargetAtom().

◆ SelectRandomAtomNumber()

G4int G4VEmModel::SelectRandomAtomNumber ( const G4Material mat)

Definition at line 315 of file G4VEmModel.cc.

316{
317 // this algorith assumes that cross section is proportional to
318 // number electrons multiplied by number of atoms
319 size_t nn = mat->GetNumberOfElements();
320 fCurrentElement = mat->GetElement(0);
321 if(1 < nn) {
322 const G4double* at = mat->GetVecNbOfAtomsPerVolume();
324 for(size_t i=0; i<nn; ++i) {
325 tot -= at[i];
326 if(tot <= 0.0) {
327 fCurrentElement = mat->GetElement(i);
328 break;
329 }
330 }
331 }
332 return fCurrentElement->GetZasInt();
333}
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207

Referenced by G4AtimaEnergyLossModel::SampleSecondaries(), G4BetheBlochModel::SampleSecondaries(), G4BraggIonModel::SampleSecondaries(), G4BraggModel::SampleSecondaries(), G4ICRU73QOModel::SampleSecondaries(), G4LindhardSorensenIonModel::SampleSecondaries(), G4MollerBhabhaModel::SampleSecondaries(), and G4IonParametrisedLossModel::SampleSecondaries().

◆ SelectTargetAtom()

◆ SetActivationHighEnergyLimit()

void G4VEmModel::SetActivationHighEnergyLimit ( G4double  val)
inline

◆ SetActivationLowEnergyLimit()

◆ SetAngularDistribution()

void G4VEmModel::SetAngularDistribution ( G4VEmAngularDistribution p)
inline

Definition at line 618 of file G4VEmModel.hh.

619{
620 if(p != anglModel) {
621 delete anglModel;
622 anglModel = p;
623 }
624}

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4BetheHeitlerModel::G4BetheHeitlerModel(), G4DNABornIonisationModel1::G4DNABornIonisationModel1(), G4DNABornIonisationModel2::G4DNABornIonisationModel2(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DNARuddIonisationModel::G4DNARuddIonisationModel(), G4eBremParametrizedModel::G4eBremParametrizedModel(), G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel(), G4LivermoreIonisationModel::G4LivermoreIonisationModel(), G4LivermorePhotoElectricModel::G4LivermorePhotoElectricModel(), G4LivermorePolarizedPhotoElectricModel::G4LivermorePolarizedPhotoElectricModel(), G4LivermoreRayleighModel::G4LivermoreRayleighModel(), G4MicroElecInelasticModel::G4MicroElecInelasticModel(), G4MicroElecInelasticModel_new::G4MicroElecInelasticModel_new(), G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(), G4MuPairProductionModel::G4MuPairProductionModel(), G4PAIModel::G4PAIModel(), G4PAIPhotModel::G4PAIPhotModel(), G4PairProductionRelModel::G4PairProductionRelModel(), G4PEEffectFluoModel::G4PEEffectFluoModel(), G4SeltzerBergerModel::G4SeltzerBergerModel(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), and G4MollerBhabhaModel::Initialise().

◆ SetAngularGeneratorFlag()

void G4VEmModel::SetAngularGeneratorFlag ( G4bool  val)
inline

Definition at line 715 of file G4VEmModel.hh.

716{
717 useAngularGenerator = val;
718}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ SetCrossSectionTable()

void G4VEmModel::SetCrossSectionTable ( G4PhysicsTable p,
G4bool  isLocal 
)

Definition at line 464 of file G4VEmModel.cc.

465{
466 if(p != xSectionTable) {
467 if(xSectionTable != nullptr && localTable) {
469 delete xSectionTable;
470 }
471 xSectionTable = p;
472 }
473 localTable = isLocal;
474}

Referenced by G4VMultipleScattering::BuildPhysicsTable().

◆ SetCurrentCouple()

◆ SetCurrentElement()

void G4VEmModel::SetCurrentElement ( const G4Element elm)
inlineprotected

Definition at line 487 of file G4VEmModel.hh.

488{
489 fCurrentElement = elm;
490 fCurrentIsotope = nullptr;
491}

Referenced by ComputeCrossSectionPerAtom(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), and SelectIsotopeNumber().

◆ SetDeexcitationFlag()

void G4VEmModel::SetDeexcitationFlag ( G4bool  val)
inline

Definition at line 813 of file G4VEmModel.hh.

814{
815 flagDeexcitation = val;
816}

Referenced by G4DNABornIonisationModel1::G4DNABornIonisationModel1(), G4DNABornIonisationModel2::G4DNABornIonisationModel2(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DNARuddIonisationModel::G4DNARuddIonisationModel(), G4KleinNishinaModel::G4KleinNishinaModel(), G4LEPTSIonisationModel::G4LEPTSIonisationModel(), G4LivermoreComptonModel::G4LivermoreComptonModel(), G4LivermoreComptonModifiedModel::G4LivermoreComptonModifiedModel(), G4LivermorePhotoElectricModel::G4LivermorePhotoElectricModel(), G4LivermorePolarizedComptonModel::G4LivermorePolarizedComptonModel(), G4LivermorePolarizedPhotoElectricGDModel::G4LivermorePolarizedPhotoElectricGDModel(), G4LivermorePolarizedPhotoElectricModel::G4LivermorePolarizedPhotoElectricModel(), G4LowEPComptonModel::G4LowEPComptonModel(), G4LowEPPolarizedComptonModel::G4LowEPPolarizedComptonModel(), G4MicroElecInelasticModel::G4MicroElecInelasticModel(), G4MicroElecInelasticModel_new::G4MicroElecInelasticModel_new(), G4PEEffectFluoModel::G4PEEffectFluoModel(), G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel(), G4PenelopeComptonModel::G4PenelopeComptonModel(), G4PenelopeIonisationModel::G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), and G4LindhardSorensenIonModel::Initialise().

◆ SetElementSelectors()

void G4VEmModel::SetElementSelectors ( std::vector< G4EmElementSelector * > *  p)
inline

◆ SetFluctuationFlag()

void G4VEmModel::SetFluctuationFlag ( G4bool  val)
inline

Definition at line 722 of file G4VEmModel.hh.

723{
724 lossFlucFlag = val;
725}

Referenced by G4EmCalculator::ComputeNuclearDEDX().

◆ SetForceBuildTable()

void G4VEmModel::SetForceBuildTable ( G4bool  val)
inline

Definition at line 820 of file G4VEmModel.hh.

821{
822 flagForceBuildTable = val;
823}

◆ SetHighEnergyLimit()

void G4VEmModel::SetHighEnergyLimit ( G4double  val)
inline

Definition at line 757 of file G4VEmModel.hh.

758{
759 highLimit = val;
760}

Referenced by LBE::ConstructEM(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmDNAPhysics_option7::ConstructProcess(), G4EmLEPTSPhysics::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4BraggIonModel::G4BraggIonModel(), G4BraggModel::G4BraggModel(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNAIonElasticModel::G4DNAIonElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel(), G4DNATransformElectronModel::G4DNATransformElectronModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel(), G4ICRU73QOModel::G4ICRU73QOModel(), G4MicroElecElasticModel::G4MicroElecElasticModel(), G4MicroElecElasticModel_new::G4MicroElecElasticModel_new(), G4PenelopeAnnihilationModel::G4PenelopeAnnihilationModel(), G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel(), G4PenelopeComptonModel::G4PenelopeComptonModel(), G4PenelopeGammaConversionModel::G4PenelopeGammaConversionModel(), G4PenelopeIonisationModel::G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel(), G4PenelopeRayleighModel::G4PenelopeRayleighModel(), G4PenelopeRayleighModelMI::G4PenelopeRayleighModelMI(), G4TDNAOneStepThermalizationModel< MODEL >::G4TDNAOneStepThermalizationModel(), G4XrayRayleighModel::G4XrayRayleighModel(), G4VLEPTSModel::Init(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNAModelInterface::Initialise(), G4DNAIonElasticModel::Initialise(), G4hBremsstrahlung::InitialiseEnergyLossProcess(), G4hhIonisation::InitialiseEnergyLossProcess(), G4hPairProduction::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4ePairProduction::InitialiseEnergyLossProcess(), G4MuBremsstrahlung::InitialiseEnergyLossProcess(), G4MuIonisation::InitialiseEnergyLossProcess(), G4MuPairProduction::InitialiseEnergyLossProcess(), G4ePolarizedBremsstrahlung::InitialiseEnergyLossProcess(), G4ePolarizedIonisation::InitialiseEnergyLossProcess(), G4alphaIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4eIonisation::InitialiseEnergyLossProcess(), G4hIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectricEffect::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit(), G4VEmAdjointModel::SetHighEnergyLimit(), G4mplIonisationModel::SetParticle(), and G4mplIonisationWithDeltaModel::SetParticle().

◆ SetLocked()

void G4VEmModel::SetLocked ( G4bool  val)
inline

Definition at line 874 of file G4VEmModel.hh.

875{
876 isLocked = val;
877}

◆ SetLowEnergyLimit()

void G4VEmModel::SetLowEnergyLimit ( G4double  val)
inline

Definition at line 764 of file G4VEmModel.hh.

765{
766 lowLimit = val;
767}

Referenced by G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4DNASancheExcitationModel::ExtendLowEnergyLimit(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4BetheBlochModel::G4BetheBlochModel(), G4BetheHeitler5DModel::G4BetheHeitler5DModel(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNAIonElasticModel::G4DNAIonElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel(), G4DNATransformElectronModel::G4DNATransformElectronModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4DummyModel::G4DummyModel(), G4eBremParametrizedModel::G4eBremParametrizedModel(), G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel(), G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel(), G4LindhardSorensenIonModel::G4LindhardSorensenIonModel(), G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel(), G4MicroElecElasticModel::G4MicroElecElasticModel(), G4MicroElecElasticModel_new::G4MicroElecElasticModel_new(), G4SeltzerBergerModel::G4SeltzerBergerModel(), G4TDNAOneStepThermalizationModel< MODEL >::G4TDNAOneStepThermalizationModel(), G4VLEPTSModel::Init(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4BetheHeitler5DModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNAModelInterface::Initialise(), G4DNAIonElasticModel::Initialise(), G4hBremsstrahlung::InitialiseEnergyLossProcess(), G4hhIonisation::InitialiseEnergyLossProcess(), G4hPairProduction::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4ePairProduction::InitialiseEnergyLossProcess(), G4MuBremsstrahlung::InitialiseEnergyLossProcess(), G4MuIonisation::InitialiseEnergyLossProcess(), G4MuPairProduction::InitialiseEnergyLossProcess(), G4ePolarizedBremsstrahlung::InitialiseEnergyLossProcess(), G4ePolarizedIonisation::InitialiseEnergyLossProcess(), G4alphaIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4eIonisation::InitialiseEnergyLossProcess(), G4hIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectricEffect::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4VEmAdjointModel::SetLowEnergyLimit(), G4mplIonisationModel::SetParticle(), and G4mplIonisationWithDeltaModel::SetParticle().

◆ SetLPMFlag()

◆ SetMasterThread()

void G4VEmModel::SetMasterThread ( G4bool  val)
inline

◆ SetParticleChange()

◆ SetPolarAngleLimit()

void G4VEmModel::SetPolarAngleLimit ( G4double  val)
inline

Definition at line 792 of file G4VEmModel.hh.

793{
794 if(!isLocked) { polarAngleLimit = val; }
795}

Referenced by G4CoulombScattering::InitialiseProcess(), G4VEmProcess::PreparePhysicsTable(), and G4VMultipleScattering::PreparePhysicsTable().

◆ SetSecondaryThreshold()

◆ SetTripletModel()

void G4VEmModel::SetTripletModel ( G4VEmModel p)
inline

Definition at line 635 of file G4VEmModel.hh.

636{
637 if(p != fTripletModel) {
638 delete fTripletModel;
639 fTripletModel = p;
640 }
641}

Referenced by G4eplusTo2GammaOKVIModel::G4eplusTo2GammaOKVIModel().

◆ SetupForMaterial()

◆ SetUseBaseMaterials()

void G4VEmModel::SetUseBaseMaterials ( G4bool  val)
inline

Definition at line 743 of file G4VEmModel.hh.

744{
745 useBaseMaterials = val;
746}

Referenced by G4VMscModel::G4VMscModel().

◆ StartTracking()

void G4VEmModel::StartTracking ( G4Track )
virtual

◆ UseAngularGeneratorFlag()

◆ UseBaseMaterials()

G4bool G4VEmModel::UseBaseMaterials ( ) const
inline

Definition at line 750 of file G4VEmModel.hh.

751{
752 return useBaseMaterials;
753}

Referenced by G4LossTableBuilder::BuildTableForModel().

◆ Value()

Member Data Documentation

◆ fElementData

G4ElementData* G4VEmModel::fElementData
protected

◆ idxTable

size_t G4VEmModel::idxTable
protected

◆ inveplus

G4double G4VEmModel::inveplus
protected

◆ lossFlucFlag

G4bool G4VEmModel::lossFlucFlag
protected

Definition at line 445 of file G4VEmModel.hh.

Referenced by SetFluctuationFlag().

◆ pBaseMaterial

const G4Material* G4VEmModel::pBaseMaterial
protected

◆ pFactor

G4double G4VEmModel::pFactor
protected

Definition at line 447 of file G4VEmModel.hh.

Referenced by ComputeDEDX(), CrossSection(), SetCurrentCouple(), and Value().

◆ pParticleChange

◆ theDensityFactor

const std::vector<G4double>* G4VEmModel::theDensityFactor
protected

Definition at line 442 of file G4VEmModel.hh.

Referenced by G4VEmModel().

◆ theDensityIdx

const std::vector<G4int>* G4VEmModel::theDensityIdx
protected

Definition at line 443 of file G4VEmModel.hh.

Referenced by G4VEmModel().

◆ xSectionTable


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