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

#include <G4VEmProcess.hh>

Inheritance diagram for G4VEmProcess:

Public Member Functions

 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 ~G4VEmProcess () override
void ProcessDescription (std::ostream &outFile) const override
void PreparePhysicsTable (const G4ParticleDefinition &) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
void StartTracking (G4Track *) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
virtual G4VEmProcessGetEmProcess (const G4String &name)
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
G4double GetCrossSection (const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4double MeanFreePath (const G4Track &track)
void SetLambdaBinning (G4int nbins)
void SetMinKinEnergy (G4double e)
void SetMinKinEnergyPrim (G4double e)
void SetMaxKinEnergy (G4double e)
G4PhysicsTableLambdaTable () const
G4PhysicsTableLambdaTablePrim () const
void SetLambdaTable (G4PhysicsTable *)
void SetLambdaTablePrim (G4PhysicsTable *)
std::vector< G4double > * EnergyOfCrossSectionMax () const
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
G4CrossSectionType CrossSectionType () const
void SetCrossSectionType (G4CrossSectionType val)
const G4ParticleDefinitionParticle () const
const G4ParticleDefinitionSecondaryParticle () const
G4VEmModelSelectModelForMaterial (G4double kinEnergy, std::size_t idxCouple) const
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel (G4VEmModel *, G4int index=0)
G4int NumberOfModels () const
G4VEmModelEmModel (std::size_t index=0) const
const G4VEmModelGetCurrentModel () const
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
const G4ElementGetCurrentElement () const
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
G4double CrossSectionBiasingFactor () const
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
void SetEmMasterProcess (const G4VEmProcess *)
void SetBuildTableFlag (G4bool val)
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
G4bool UseBaseMaterial () const
void BuildLambdaTable ()
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 G4VEmProcess (G4VEmProcess &)=delete
G4VEmProcessoperator= (const G4VEmProcess &right)=delete
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 G4VProcess (const G4VProcess &right)
virtual ~G4VProcess ()
G4VProcessoperator= (const G4VProcess &)=delete
G4bool operator== (const G4VProcess &right) const
G4bool operator!= (const G4VProcess &right) const
G4double GetCurrentInteractionLength () const
void SetPILfactor (G4double value)
G4double GetPILfactor () const
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4bool IsApplicable (const G4ParticleDefinition &)
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
const G4StringGetProcessName () const
G4ProcessType GetProcessType () const
void SetProcessType (G4ProcessType)
G4int GetProcessSubType () const
void SetProcessSubType (G4int)
virtual const G4VProcessGetCreatorProcess () const
virtual void EndTracking ()
virtual void SetProcessManager (const G4ProcessManager *)
virtual const G4ProcessManagerGetProcessManager ()
virtual void ResetNumberOfInteractionLengthLeft ()
G4double GetNumberOfInteractionLengthLeft () const
G4double GetTotalNumberOfInteractionLengthTraversed () const
G4bool isAtRestDoItIsEnabled () const
G4bool isAlongStepDoItIsEnabled () const
G4bool isPostStepDoItIsEnabled () const
virtual void DumpInfo () const
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const
virtual void SetMasterProcess (G4VProcess *masterP)
const G4VProcessGetMasterProcess () const
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)

Protected Member Functions

virtual void StreamProcessInfo (std::ostream &) const
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
G4VEmModelSelectModel (G4double kinEnergy, std::size_t)
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
void DefineMaterial (const G4MaterialCutsCouple *couple)
G4int LambdaBinning () const
G4double MinKinEnergy () const
G4double MaxKinEnergy () const
G4double PolarAngleLimit () const
G4ParticleChangeForGammaGetParticleChange ()
void SetParticle (const G4ParticleDefinition *p)
void SetSecondaryParticle (const G4ParticleDefinition *p)
std::size_t CurrentMaterialCutsCoupleIndex () const
const G4MaterialCutsCoupleMaterialCutsCouple () const
G4bool ApplyCuts () const
G4double GetGammaEnergyCut ()
G4double GetElectronEnergyCut ()
void SetStartFromNullFlag (G4bool val)
void SetSplineFlag (G4bool val)
const G4ElementGetTargetElement () const
const G4IsotopeGetTargetIsotope () const
G4int DensityIndex (G4int idx) const
G4double DensityFactor (G4int idx) const
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

const G4MaterialCutsCouplecurrentCouple = nullptr
const G4MaterialcurrentMaterial = nullptr
G4EmBiasingManagerbiasManager = nullptr
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
G4double mfpKinEnergy = DBL_MAX
G4double preStepKinEnergy = 0.0
G4double preStepLambda = 0.0
G4int mainSecondaries = 1
G4int secID = _EM
G4int fluoID = _Fluorescence
G4int augerID = _AugerElectron
G4int biasID = _EM
G4int tripletID = _TripletElectron
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
std::size_t coupleIdxLambda = 0
std::size_t idxLambda = 0
G4bool isTheMaster = false
G4bool baseMat = false
std::vector< G4DynamicParticle * > secParticles
G4ParticleChangeForGamma fParticleChange
Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
G4VParticleChangepParticleChange = nullptr
G4ParticleChange aParticleChange
G4double theNumberOfInteractionLengthLeft = -1.0
G4double currentInteractionLength = -1.0
G4double theInitialNumberOfInteractionLength = -1.0
G4String theProcessName
G4String thePhysicsTableFileName
G4ProcessType theProcessType = fNotDefined
G4int theProcessSubType = -1
G4double thePILfactor = 1.0
G4int verboseLevel = 0
G4bool enableAtRestDoIt = true
G4bool enableAlongStepDoIt = true
G4bool enablePostStepDoIt = true

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)

Detailed Description

Definition at line 76 of file G4VEmProcess.hh.

Constructor & Destructor Documentation

◆ G4VEmProcess() [1/2]

G4VEmProcess::G4VEmProcess ( const G4String & name,
G4ProcessType type = fElectromagnetic )

Definition at line 79 of file G4VEmProcess.cc.

79 :
80 G4VDiscreteProcess(name, type)
81{
82 theParameters = G4EmParameters::Instance();
84
85 // Size of tables
86 minKinEnergy = 0.1*CLHEP::keV;
87 maxKinEnergy = 100.0*CLHEP::TeV;
88
89 // default lambda factor
90 invLambdaFactor = 1.0/lambdaFactor;
91
92 // particle types
93 theGamma = G4Gamma::Gamma();
94 theElectron = G4Electron::Electron();
95 thePositron = G4Positron::Positron();
96
98 fParticleChange.SetSecondaryWeightByProcess(true);
99 secParticles.reserve(5);
100
101 modelManager = new G4EmModelManager();
102 lManager = G4LossTableManager::Instance();
103 lManager->Register(this);
104 isTheMaster = lManager->IsMaster();
105 G4LossTableBuilder* bld = lManager->GetTableBuilder();
106 theDensityFactor = bld->GetDensityFactors();
107 theDensityIdx = bld->GetCoupleIndexes();
108}
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static const std::vector< G4double > * GetDensityFactors()
static const std::vector< G4int > * GetCoupleIndexes()
static G4LossTableManager * Instance()
static G4Positron * Positron()
Definition G4Positron.cc:90
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
std::vector< G4DynamicParticle * > secParticles
G4ParticleChangeForGamma fParticleChange
void SetVerboseLevel(G4int value)
G4VParticleChange * pParticleChange

Referenced by G4GammaGeneralProcess::AddEmProcess(), BuildPhysicsTable(), G4ComptonScattering::G4ComptonScattering(), G4CoulombScattering::G4CoulombScattering(), G4DNAAttachment::G4DNAAttachment(), G4DNAChargeDecrease::G4DNAChargeDecrease(), G4DNAChargeIncrease::G4DNAChargeIncrease(), G4DNADissociation::G4DNADissociation(), G4DNADoubleIonisation::G4DNADoubleIonisation(), G4DNAElastic::G4DNAElastic(), G4DNAElectronSolvation::G4DNAElectronSolvation(), G4DNAExcitation::G4DNAExcitation(), G4DNAIonisation::G4DNAIonisation(), G4DNAPlasmonExcitation::G4DNAPlasmonExcitation(), G4DNAPositronium::G4DNAPositronium(), G4DNAQuadrupleIonisation::G4DNAQuadrupleIonisation(), G4DNARotExcitation::G4DNARotExcitation(), G4DNATripleIonisation::G4DNATripleIonisation(), G4DNAVibExcitation::G4DNAVibExcitation(), G4eeToHadrons::G4eeToHadrons(), G4eplusAnnihilation::G4eplusAnnihilation(), G4GammaConversion::G4GammaConversion(), G4GammaGeneralProcess::G4GammaGeneralProcess(), G4JAEAElasticScattering::G4JAEAElasticScattering(), G4MicroElecElastic::G4MicroElecElastic(), G4MicroElecInelastic::G4MicroElecInelastic(), G4MicroElecLOPhononScattering::G4MicroElecLOPhononScattering(), G4NuclearStopping::G4NuclearStopping(), G4PhotoElectricEffect::G4PhotoElectricEffect(), G4PolarizedCompton::G4PolarizedCompton(), G4PolarizedGammaConversion::G4PolarizedGammaConversion(), G4PolarizedPhotoElectric::G4PolarizedPhotoElectric(), G4RayleighScattering::G4RayleighScattering(), G4VEmProcess(), G4GammaGeneralProcess::GetEmProcess(), GetEmProcess(), InitialiseProcess(), G4GammaGeneralProcess::operator=(), operator=(), G4GammaGeneralProcess::SelectEmProcess(), and SetEmMasterProcess().

◆ ~G4VEmProcess()

G4VEmProcess::~G4VEmProcess ( )
override

Definition at line 112 of file G4VEmProcess.cc.

113{
114 if (isTheMaster) {
116 }
117 delete modelManager;
118 delete biasManager;
119 lManager->DeRegister(this);
120}
G4EmBiasingManager * biasManager
std::vector< G4double > * theEnergyOfCrossSectionMax

◆ G4VEmProcess() [2/2]

G4VEmProcess::G4VEmProcess ( G4VEmProcess & )
delete

Member Function Documentation

◆ ActivateForcedInteraction()

void G4VEmProcess::ActivateForcedInteraction ( G4double length = 0.0,
const G4String & r = "",
G4bool flag = true )

Definition at line 776 of file G4VEmProcess.cc.

778{
779 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
780 if(1 < verboseLevel) {
781 G4cout << "### ActivateForcedInteraction: for "
782 << particle->GetParticleName()
783 << " and process " << GetProcessName()
784 << " length(mm)= " << length/mm
785 << " in G4Region <" << r
786 << "> weightFlag= " << flag
787 << G4endl;
788 }
789 weightFlag = flag;
790 biasManager->ActivateForcedInteraction(length, r);
791}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
const G4String & GetProcessName() const

Referenced by G4EmExtraParameters::DefineRegParamForEM().

◆ ActivateSecondaryBiasing()

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

Definition at line 796 of file G4VEmProcess.cc.

799{
800 if (0.0 <= factor) {
801
802 // Range cut can be applied only for e-
803 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
804 { return; }
805
806 if(!biasManager) { biasManager = new G4EmBiasingManager(); }
807 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
808 if(1 < verboseLevel) {
809 G4cout << "### ActivateSecondaryBiasing: for "
810 << " process " << GetProcessName()
811 << " factor= " << factor
812 << " in G4Region <" << region
813 << "> energyLimit(MeV)= " << energyLimit/MeV
814 << G4endl;
815 }
816 }
817}

Referenced by G4EmExtraParameters::DefineRegParamForEM().

◆ AddEmModel()

void G4VEmProcess::AddEmModel ( G4int order,
G4VEmModel * ptr,
const G4Region * region = nullptr )

Definition at line 124 of file G4VEmProcess.cc.

126{
127 if(nullptr == ptr) { return; }
128 G4VEmFluctuationModel* fm = nullptr;
129 modelManager->AddEmModel(order, ptr, fm, region);
131}
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)

Referenced by G4EmDNABuilder::ConstructDNAElectronPhysics(), LBE::ConstructEM(), G4EmLivermorePhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNADoubleIonisation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAElectronSolvation::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPlasmonExcitation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNAQuadrupleIonisation::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNATripleIonisation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4eeToHadrons::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4JAEAElasticScattering::InitialiseProcess(), G4MicroElecElastic::InitialiseProcess(), G4MicroElecInelastic::InitialiseProcess(), G4MicroElecLOPhononScattering::InitialiseProcess(), G4NuclearStopping::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4RayleighScattering::InitialiseProcess(), and G4EmConfigurator::PrepareModels().

◆ ApplyCuts()

G4bool G4VEmProcess::ApplyCuts ( ) const
inlineprotected

Definition at line 613 of file G4VEmProcess.hh.

614{
615 return applyCuts;
616}

◆ BuildLambdaTable()

void G4VEmProcess::BuildLambdaTable ( )

Definition at line 245 of file G4VEmProcess.cc.

246{
247 G4double scale = theParameters->MaxKinEnergy()/theParameters->MinKinEnergy();
248 G4int nbin =
249 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
250 if(actBinning) { nbin = std::max(nbin, nLambdaBins); }
251 scale = nbin/G4Log(scale);
252
253 G4LossTableBuilder* bld = lManager->GetTableBuilder();
254 G4EmTableUtil::BuildLambdaTable(this, particle, modelManager,
255 bld, theLambdaTable, theLambdaTablePrim,
256 minKinEnergy, minKinEnergyPrim,
257 maxKinEnergy, scale, verboseLevel,
258 startFromNull, splineFlag);
259}
G4double G4Log(G4double x)
Definition G4Log.hh:169
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static void BuildLambdaTable(G4VEmProcess *proc, const G4ParticleDefinition *part, G4EmModelManager *modelManager, G4LossTableBuilder *bld, G4PhysicsTable *theLambdaTable, G4PhysicsTable *theLambdaTablePrim, const G4double minKinEnergy, const G4double minKinEnergyPrim, const G4double maxKinEnergy, const G4double scale, const G4int verbose, const G4bool startFromNull, const G4bool splineFlag)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4EmTableUtil::BuildEmProcess().

◆ BuildPhysicsTable()

void G4VEmProcess::BuildPhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 228 of file G4VEmProcess.cc.

229{
230 if(nullptr == masterProc) {
231 if(isTheMaster) { masterProc = this; }
232 else { masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());}
233 }
234 G4int nModels = modelManager->NumberOfModels();
235 G4bool isLocked = theParameters->IsPrintLocked();
236 G4bool toBuild = (buildLambdaTable || minKinEnergyPrim < maxKinEnergy);
237
238 G4EmTableUtil::BuildEmProcess(this, masterProc, particle, &part,
239 nModels, verboseLevel, isTheMaster,
240 isLocked, toBuild, baseMat);
241}
bool G4bool
Definition G4Types.hh:86
static void BuildEmProcess(G4VEmProcess *proc, const G4VEmProcess *masterProc, const G4ParticleDefinition *firstPart, const G4ParticleDefinition *part, const G4int nModels, const G4int verb, const G4bool master, const G4bool isLocked, const G4bool toBuild, G4bool &baseMat)
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
const G4VProcess * GetMasterProcess() const

Referenced by G4PolarizedAnnihilation::BuildPhysicsTable(), and G4PolarizedCompton::BuildPhysicsTable().

◆ ComputeCrossSectionPerAtom()

G4double G4VEmProcess::ComputeCrossSectionPerAtom ( G4double kineticEnergy,
G4double Z,
G4double A = 0.,
G4double cut = 0.0 )

Definition at line 712 of file G4VEmProcess.cc.

714{
716 return (currentModel) ?
717 currentModel->ComputeCrossSectionPerAtom(currentParticle, kinEnergy,
718 Z, A, cut) : 0.0;
719}
const G4double A[17]
std::size_t currentCoupleIndex
G4VEmModel * SelectModel(G4double kinEnergy, std::size_t)

◆ CrossSectionBiasingFactor()

G4double G4VEmProcess::CrossSectionBiasingFactor ( ) const
inline

Definition at line 641 of file G4VEmProcess.hh.

642{
643 return biasFactor;
644}

◆ CrossSectionType()

G4CrossSectionType G4VEmProcess::CrossSectionType ( ) const
inline

Definition at line 712 of file G4VEmProcess.hh.

713{
714 return fXSType;
715}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ CurrentMaterialCutsCoupleIndex()

std::size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex ( ) const
inlineprotected

Definition at line 443 of file G4VEmProcess.hh.

444{
445 return currentCoupleIndex;
446}

◆ CurrentSetup()

void G4VEmProcess::CurrentSetup ( const G4MaterialCutsCouple * couple,
G4double energy )
inline

Definition at line 584 of file G4VEmProcess.hh.

585{
586 DefineMaterial(couple);
587 SelectModel(energy*massRatio, currentCoupleIndex);
588}
void DefineMaterial(const G4MaterialCutsCouple *couple)

Referenced by GetCrossSection(), GetLambda(), MeanFreePath(), and G4GammaGeneralProcess::SelectEmProcess().

◆ DefineMaterial()

void G4VEmProcess::DefineMaterial ( const G4MaterialCutsCouple * couple)
inlineprotected

Definition at line 471 of file G4VEmProcess.hh.

472{
473 if (couple != currentCouple) {
474 currentCouple = couple;
475 baseMaterial = currentMaterial = couple->GetMaterial();
477 fFactor = biasFactor;
479 if (baseMat) {
480 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
481 if (nullptr != currentMaterial->GetBaseMaterial())
482 baseMaterial = currentMaterial->GetBaseMaterial();
483 fFactor *= (*theDensityFactor)[currentCoupleIndex];
484 }
485 }
486}
const G4Material * GetMaterial() const
G4double mfpKinEnergy
std::size_t basedCoupleIndex
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial
#define DBL_MAX
Definition templates.hh:62

Referenced by G4eplusAnnihilation::AtRestDoIt(), CurrentSetup(), LambdaPhysicsVector(), and PostStepGetPhysicalInteractionLength().

◆ DensityFactor()

G4double G4VEmProcess::DensityFactor ( G4int idx) const
inlineprotected

Definition at line 769 of file G4VEmProcess.hh.

770{
771 return (*theDensityFactor)[idx];
772}

Referenced by G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength().

◆ DensityIndex()

G4int G4VEmProcess::DensityIndex ( G4int idx) const
inlineprotected

Definition at line 762 of file G4VEmProcess.hh.

763{
764 return (*theDensityIdx)[idx];
765}

Referenced by G4GammaGeneralProcess::BuildPhysicsTable(), and G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength().

◆ EmModel()

G4VEmModel * G4VEmProcess::EmModel ( std::size_t index = 0) const
inline

Definition at line 804 of file G4VEmProcess.hh.

805{
806 return (index < emModels.size()) ? emModels[index] : nullptr;
807}

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4EmCalculator::ComputeNuclearDEDX(), G4EmDNAChemistry::ConstructProcess(), G4EmDNAChemistry_option1::ConstructProcess(), G4EmDNAChemistry_option2::ConstructProcess(), G4EmDNAChemistry_option3::ConstructProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNADoubleIonisation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAElectronSolvation::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPlasmonExcitation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNAQuadrupleIonisation::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNATripleIonisation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4JAEAElasticScattering::InitialiseProcess(), G4MicroElecElastic::InitialiseProcess(), G4MicroElecInelastic::InitialiseProcess(), G4MicroElecLOPhononScattering::InitialiseProcess(), G4NuclearStopping::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4RayleighScattering::InitialiseProcess(), G4DNAAttachment::PrintInfo(), G4DNADissociation::PrintInfo(), G4DNADoubleIonisation::PrintInfo(), G4DNAElastic::PrintInfo(), G4DNAIonisation::PrintInfo(), G4DNAPlasmonExcitation::PrintInfo(), G4DNAPositronium::PrintInfo(), G4DNAQuadrupleIonisation::PrintInfo(), G4DNATripleIonisation::PrintInfo(), G4DNAChargeDecrease::ProcessDescription(), and G4DNAChargeIncrease::ProcessDescription().

◆ EnergyOfCrossSectionMax()

std::vector< G4double > * G4VEmProcess::EnergyOfCrossSectionMax ( ) const
inline

Definition at line 676 of file G4VEmProcess.hh.

677{
679}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ GetCrossSection()

G4double G4VEmProcess::GetCrossSection ( const G4double kinEnergy,
const G4MaterialCutsCouple * couple )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 692 of file G4VEmProcess.cc.

694{
695 CurrentSetup(couple, kinEnergy);
696 return GetCurrentLambda(kinEnergy, G4Log(kinEnergy));
697}
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)

Referenced by G4EmCalculator::GetCrossSectionPerVolume().

◆ GetCurrentElement()

const G4Element * G4VEmProcess::GetCurrentElement ( ) const

Definition at line 734 of file G4VEmProcess.cc.

735{
736 return (nullptr != currentModel) ?
737 currentModel->GetCurrentElement(currentMaterial) : nullptr;
738}

Referenced by GetTargetIsotope().

◆ GetCurrentModel()

const G4VEmModel * G4VEmProcess::GetCurrentModel ( ) const
inline

Definition at line 783 of file G4VEmProcess.hh.

784{
785 return currentModel;
786}

◆ GetElectronEnergyCut()

G4double G4VEmProcess::GetElectronEnergyCut ( )
inlineprotected

Definition at line 464 of file G4VEmProcess.hh.

465{
466 return (*theCutsElectron)[currentCoupleIndex];
467}

◆ GetEmProcess()

G4VEmProcess * G4VEmProcess::GetEmProcess ( const G4String & name)
virtual

Reimplemented in G4GammaGeneralProcess.

Definition at line 867 of file G4VEmProcess.cc.

868{
869 return (nam == GetProcessName()) ? this : nullptr;
870}

◆ GetGammaEnergyCut()

G4double G4VEmProcess::GetGammaEnergyCut ( )
inlineprotected

Definition at line 457 of file G4VEmProcess.hh.

458{
459 return (*theCutsGamma)[currentCoupleIndex];
460}

Referenced by G4eplusAnnihilation::AtRestDoIt().

◆ GetLambda()

G4double G4VEmProcess::GetLambda ( G4double kinEnergy,
const G4MaterialCutsCouple * couple,
G4double logKinEnergy )
inline

Definition at line 593 of file G4VEmProcess.hh.

595{
596 CurrentSetup(couple, kinEnergy);
597 return GetCurrentLambda(kinEnergy, logKinEnergy);
598}

◆ GetMeanFreePath()

G4double G4VEmProcess::GetMeanFreePath ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VDiscreteProcess.

Definition at line 701 of file G4VEmProcess.cc.

704{
706 return G4VEmProcess::MeanFreePath(track);
707}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double MeanFreePath(const G4Track &track)

Referenced by G4PolarizedAnnihilation::GetMeanFreePath(), and G4PolarizedCompton::GetMeanFreePath().

◆ GetModelByIndex()

G4VEmModel * G4VEmProcess::GetModelByIndex ( G4int idx = 0,
G4bool ver = false ) const
inline

Definition at line 811 of file G4VEmProcess.hh.

812{
813 return modelManager->GetModel(idx, ver);
814}

Referenced by G4EmTableUtil::BuildEmProcess(), and StartTracking().

◆ GetParticleChange()

G4ParticleChangeForGamma * G4VEmProcess::GetParticleChange ( )
inlineprotected

Definition at line 726 of file G4VEmProcess.hh.

727{
728 return &fParticleChange;
729}

◆ GetTargetElement()

const G4Element * G4VEmProcess::GetTargetElement ( ) const
protected

Definition at line 742 of file G4VEmProcess.cc.

743{
744 return (nullptr != currentModel) ?
745 currentModel->GetCurrentElement(currentMaterial) : nullptr;
746}

◆ GetTargetIsotope()

const G4Isotope * G4VEmProcess::GetTargetIsotope ( ) const
protected

Definition at line 750 of file G4VEmProcess.cc.

751{
752 return (nullptr != currentModel) ?
753 currentModel->GetCurrentIsotope(GetCurrentElement()) : nullptr;
754}
const G4Element * GetCurrentElement() const

◆ InitialiseProcess()

◆ LambdaBinning()

G4int G4VEmProcess::LambdaBinning ( ) const
inlineprotected

Definition at line 620 of file G4VEmProcess.hh.

621{
622 return nLambdaBins;
623}

◆ LambdaPhysicsVector()

G4PhysicsVector * G4VEmProcess::LambdaPhysicsVector ( const G4MaterialCutsCouple * couple)
protected

Definition at line 724 of file G4VEmProcess.cc.

725{
726 DefineMaterial(couple);
727 G4PhysicsVector* newv = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy,
728 nLambdaBins, splineFlag);
729 return newv;
730}

◆ LambdaTable()

G4PhysicsTable * G4VEmProcess::LambdaTable ( ) const
inline

Definition at line 648 of file G4VEmProcess.hh.

649{
650 return theLambdaTable;
651}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ LambdaTablePrim()

G4PhysicsTable * G4VEmProcess::LambdaTablePrim ( ) const
inline

Definition at line 655 of file G4VEmProcess.hh.

656{
657 return theLambdaTablePrim;
658}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ MaterialCutsCouple()

const G4MaterialCutsCouple * G4VEmProcess::MaterialCutsCouple ( ) const
inlineprotected

Definition at line 450 of file G4VEmProcess.hh.

451{
452 return currentCouple;
453}

◆ MaxKinEnergy()

G4double G4VEmProcess::MaxKinEnergy ( ) const
inlineprotected

Definition at line 634 of file G4VEmProcess.hh.

635{
636 return maxKinEnergy;
637}

Referenced by G4eplusAnnihilation::InitialiseProcess(), and SetMinKinEnergyPrim().

◆ MeanFreePath()

G4double G4VEmProcess::MeanFreePath ( const G4Track & track)
inline

Definition at line 602 of file G4VEmProcess.hh.

603{
604 const G4double kinEnergy = track.GetKineticEnergy();
605 CurrentSetup(track.GetMaterialCutsCouple(), kinEnergy);
606 const G4double xs = GetCurrentLambda(kinEnergy,
608 return (0.0 < xs) ? 1.0/xs : DBL_MAX;
609}
G4double GetLogKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

Referenced by G4GammaGeneralProcess::GetMeanFreePath(), and GetMeanFreePath().

◆ MinKinEnergy()

G4double G4VEmProcess::MinKinEnergy ( ) const
inlineprotected

Definition at line 627 of file G4VEmProcess.hh.

628{
629 return minKinEnergy;
630}

Referenced by G4eplusAnnihilation::InitialiseProcess().

◆ NumberOfModels()

G4int G4VEmProcess::NumberOfModels ( ) const
inline

Definition at line 797 of file G4VEmProcess.hh.

798{
799 return numberOfModels;
800}

◆ operator=()

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

◆ Particle()

const G4ParticleDefinition * G4VEmProcess::Particle ( ) const
inline

Definition at line 691 of file G4VEmProcess.hh.

692{
693 return particle;
694}

◆ PolarAngleLimit()

G4double G4VEmProcess::PolarAngleLimit ( ) const
inlineprotected

Definition at line 874 of file G4VEmProcess.cc.

875{
876 return theParameters->MscThetaLimit();
877}

◆ PostStepDoIt()

G4VParticleChange * G4VEmProcess::PostStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 472 of file G4VEmProcess.cc.

474{
475 // clear number of interaction lengths in any case
478
479 fParticleChange.InitializeForPostStep(track);
480
481 // Do not make anything if particle is stopped, the annihilation then
482 // should be performed by the AtRestDoIt!
483 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
484
485 const G4double finalT = track.GetKineticEnergy();
486
487 // forced process - should happen only once per track
488 if(biasFlag) {
489 if(biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
490 biasFlag = false;
491 }
492 }
493
494 // check active and select model
495 const G4double scaledEnergy = finalT*massRatio;
496 SelectModel(scaledEnergy, currentCoupleIndex);
497 if(!currentModel->IsActive(scaledEnergy)) { return &fParticleChange; }
498
499 // Integral approach
500 if (fXSType != fEmNoIntegral) {
501 const G4double logFinalT =
503 const G4double lx = std::max(GetCurrentLambda(finalT, logFinalT), 0.0);
504#ifdef G4VERBOSE
505 if(preStepLambda < lx && 1 < verboseLevel) {
506 G4cout << "WARNING: for " << currentParticle->GetParticleName()
507 << " and " << GetProcessName() << " E(MeV)= " << finalT/MeV
508 << " preLambda= " << preStepLambda
509 << " < " << lx << " (postLambda) " << G4endl;
510 }
511#endif
512 // if false interaction then use new cross section value
513 // if both values are zero - no interaction
514 if(preStepLambda*G4UniformRand() >= lx) {
515 return &fParticleChange;
516 }
517 }
518
519 // define new weight for primary and secondaries
520 G4double weight = fParticleChange.GetParentWeight();
521 if(weightFlag) {
522 weight /= biasFactor;
523 fParticleChange.ProposeWeight(weight);
524 }
525
526#ifdef G4VERBOSE
527 if(1 < verboseLevel) {
528 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
529 << finalT/MeV
530 << " MeV; model= (" << currentModel->LowEnergyLimit()
531 << ", " << currentModel->HighEnergyLimit() << ")"
532 << G4endl;
533 }
534#endif
535
536 // sample secondaries
537 secParticles.clear();
538 currentModel->SampleSecondaries(&secParticles,
540 track.GetDynamicParticle(),
541 (*theCuts)[currentCoupleIndex]);
542
543 G4int num0 = (G4int)secParticles.size();
544
545 // splitting or Russian roulette
546 if(biasManager) {
547 if(biasManager->SecondaryBiasingRegion((G4int)currentCoupleIndex)) {
548 G4double eloss = 0.0;
549 weight *= biasManager->ApplySecondaryBiasing(
550 secParticles, track, currentModel, &fParticleChange, eloss,
552 step.GetPostStepPoint()->GetSafety());
553 if(eloss > 0.0) {
554 eloss += fParticleChange.GetLocalEnergyDeposit();
555 fParticleChange.ProposeLocalEnergyDeposit(eloss);
556 }
557 }
558 }
559
560 // save secondaries
561 G4int num = (G4int)secParticles.size();
562 if(num > 0) {
563
564 fParticleChange.SetNumberOfSecondaries(num);
565 G4double edep = fParticleChange.GetLocalEnergyDeposit();
566 G4double time = track.GetGlobalTime();
567
568 G4int n1(0), n2(0);
569 if(num0 > mainSecondaries) {
570 currentModel->FillNumberOfSecondaries(n1, n2);
571 }
572
573 for (G4int i=0; i<num; ++i) {
574 G4DynamicParticle* dp = secParticles[i];
575 if (nullptr != dp) {
576 const G4ParticleDefinition* p = dp->GetParticleDefinition();
577 G4double e = dp->GetKineticEnergy();
578 G4bool good = true;
579 if(applyCuts) {
580 if (p == theGamma) {
581 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
582
583 } else if (p == theElectron) {
584 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
585
586 } else if (p == thePositron) {
587 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
588 e < (*theCutsPositron)[currentCoupleIndex]) {
589 good = false;
590 e += 2.0*electron_mass_c2;
591 }
592 }
593 // added secondary if it is good
594 }
595 if (good) {
596 G4Track* t = new G4Track(dp, time, track.GetPosition());
598 if (biasManager) {
599 t->SetWeight(weight * biasManager->GetWeight(i));
600 } else {
601 t->SetWeight(weight);
602 }
603 pParticleChange->AddSecondary(t);
604
605 // define type of secondary
606 if(i < mainSecondaries) {
608 if(GetProcessSubType() == fComptonScattering && p == theGamma) {
610 }
611 } else if(i < mainSecondaries + n1) {
613 } else if(i < mainSecondaries + n1 + n2) {
615 } else {
616 if(i < num0) {
617 if(p == theGamma) {
619 } else {
621 }
622 } else {
624 }
625 }
626 /*
627 G4cout << "Secondary(post step) has weight " << t->GetWeight()
628 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
629 << GetProcessName() << " fluoID= " << fluoID
630 << " augerID= " << augerID <<G4endl;
631 */
632 } else {
633 delete dp;
634 edep += e;
635 }
636 }
637 }
638 fParticleChange.ProposeLocalEnergyDeposit(edep);
639 }
640
641 if(0.0 == fParticleChange.GetProposedKineticEnergy() &&
642 fAlive == fParticleChange.GetTrackStatus()) {
643 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
644 { fParticleChange.ProposeTrackStatus(fStopButAlive); }
645 else { fParticleChange.ProposeTrackStatus(fStopAndKill); }
646 }
647
648 return &fParticleChange;
649}
@ fComptonScattering
@ fEmNoIntegral
@ fAlive
@ fStopAndKill
@ fStopButAlive
#define G4UniformRand()
Definition Randomize.hh:52
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetSafety() const
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
G4int mainSecondaries
G4double preStepLambda
G4double theNumberOfInteractionLengthLeft
G4int GetProcessSubType() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4VEmProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 363 of file G4VEmProcess.cc.

367{
369 G4double x = DBL_MAX;
370
373 const G4double scaledEnergy = preStepKinEnergy*massRatio;
374 SelectModel(scaledEnergy, currentCoupleIndex);
375
376 // In models applied to ions the dynamic charge is needed
377 if (isIon) { currentModel->ChargeSquareRatio(track); }
378 /*
379 G4cout << "PostStepGetPhysicalInteractionLength: idx= " << currentCoupleIndex
380 << " couple: " << currentCouple << G4endl;
381 */
382 if(!currentModel->IsActive(scaledEnergy)) {
386 preStepLambda = 0.0;
387 return x;
388 }
389
390 // forced biasing only for primary particles
391 if(biasManager) {
392 if(0 == track.GetParentID()) {
393 if(biasFlag &&
394 biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
395 return biasManager->GetStepLimit((G4int)currentCoupleIndex, previousStepSize);
396 }
397 }
398 }
399
400 // compute mean free path
401
402 ComputeIntegralLambda(preStepKinEnergy, track);
403
404 // zero cross section
405 if(preStepLambda <= 0.0) {
408
409 } else {
410
411 // non-zero cross section
413
414 // beggining of tracking (or just after DoIt of this process)
417
418 } else {
419
421 previousStepSize/currentInteractionLength;
424 }
425
426 // new mean free path and step limit for the next step
429 }
430 return x;
431}
G4int GetParentID() const
G4double preStepKinEnergy
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength

Referenced by G4PolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), and G4PolarizedCompton::PostStepGetPhysicalInteractionLength().

◆ PreparePhysicsTable()

void G4VEmProcess::PreparePhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 146 of file G4VEmProcess.cc.

147{
148 if(nullptr == particle) { SetParticle(&part); }
149
150 if(part.GetParticleType() == "nucleus" &&
151 part.GetParticleSubType() == "generic") {
152
153 G4String pname = part.GetParticleName();
154 if(pname != "deuteron" && pname != "triton" &&
155 pname != "He3" && pname != "alpha" && pname != "alpha+" &&
156 pname != "helium" && pname != "hydrogen") {
157
158 particle = G4GenericIon::GenericIon();
159 isIon = true;
160 }
161 }
162 if(particle != &part) { return; }
163
164 lManager->PreparePhysicsTable(&part, this);
165
166 // for new run
167 currentCouple = nullptr;
168 preStepLambda = 0.0;
169 fLambdaEnergy = 0.0;
170
171 InitialiseProcess(particle);
172
173 G4LossTableBuilder* bld = lManager->GetTableBuilder();
174 const G4ProductionCutsTable* theCoupleTable=
176 theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
177 theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
178 theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
179
180 // initialisation of the process
181 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
182 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
183
184 applyCuts = theParameters->ApplyCuts();
185 lambdaFactor = theParameters->LambdaFactor();
186 invLambdaFactor = 1.0/lambdaFactor;
187 theParameters->DefineRegParamForEM(this);
188
189 // integral option may be disabled
190 if(!theParameters->Integral()) { fXSType = fEmNoIntegral; }
191
192 // prepare tables
193 if(isTheMaster) {
194 if(buildLambdaTable) {
195 if (nullptr == theData) { theData = new G4EmDataHandler(2, particle->GetParticleName()); }
196 theLambdaTable = theData->MakeTable(theLambdaTable, 0);
197 bld->InitialiseBaseMaterials(theLambdaTable);
198 }
199 // high energy table
200 if(minKinEnergyPrim < maxKinEnergy) {
201 if (nullptr == theData) { theData = new G4EmDataHandler(2, particle->GetParticleName()); }
202 theLambdaTablePrim = theData->MakeTable(theLambdaTablePrim, 1);
203 bld->InitialiseBaseMaterials(theLambdaTablePrim);
204 }
205 }
206 // models
208 numberOfModels = modelManager->NumberOfModels();
209 currentModel = modelManager->GetModel(0);
210 if(nullptr != lManager->AtomDeexcitation()) {
211 modelManager->SetFluoFlag(true);
212 }
213 // forced biasing
214 if(nullptr != biasManager) {
215 biasManager->Initialise(part, GetProcessName(), verboseLevel);
216 biasFlag = false;
217 }
218
219 theCuts =
220 G4EmTableUtil::PrepareEmProcess(this, particle, secondaryParticle,
221 modelManager, maxKinEnergy,
224}
@ idxG4ElectronCut
@ idxG4GammaCut
@ idxG4PositronCut
static const G4DataVector * PrepareEmProcess(G4VEmProcess *proc, const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4EmModelManager *modelManager, const G4double &maxKinEnergy, G4int &secID, G4int &tripletID, G4int &mainSec, const G4int &verb, const G4bool &master)
static G4GenericIon * GenericIon()
static G4bool GetBaseMaterialFlag()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
void SetParticle(const G4ParticleDefinition *p)

◆ ProcessDescription()

◆ RetrievePhysicsTable()

G4bool G4VEmProcess::RetrievePhysicsTable ( const G4ParticleDefinition * part,
const G4String & directory,
G4bool ascii )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 671 of file G4VEmProcess.cc.

674{
675 if(!isTheMaster || part != particle) { return true; }
676 G4bool yes = true;
677 if(buildLambdaTable) {
678 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTable, dir,
679 "Lambda", verboseLevel,
680 ascii, splineFlag);
681 }
682 if(yes && minKinEnergyPrim < maxKinEnergy) {
683 yes = G4EmTableUtil::RetrieveTable(this, part, theLambdaTablePrim, dir,
684 "LambdaPrim", verboseLevel,
685 ascii, splineFlag);
686 }
687 return yes;
688}
static G4bool RetrieveTable(G4VProcess *ptr, const G4ParticleDefinition *part, G4PhysicsTable *aTable, const G4String &dir, const G4String &tname, const G4int verb, const G4bool ascii, const G4bool spline)

◆ SecondaryParticle()

const G4ParticleDefinition * G4VEmProcess::SecondaryParticle ( ) const
inline

Definition at line 698 of file G4VEmProcess.hh.

699{
700 return secondaryParticle;
701}

◆ SelectModel()

G4VEmModel * G4VEmProcess::SelectModel ( G4double kinEnergy,
std::size_t  )
inlineprotected

Definition at line 491 of file G4VEmProcess.hh.

492{
493 if(1 < numberOfModels) {
494 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
495 }
496 currentModel->SetCurrentCouple(currentCouple);
497 return currentModel;
498}

Referenced by G4NuclearStopping::AlongStepDoIt(), ComputeCrossSectionPerAtom(), CurrentSetup(), PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ SelectModelForMaterial()

G4VEmModel * G4VEmProcess::SelectModelForMaterial ( G4double kinEnergy,
std::size_t idxCouple ) const
inline

Definition at line 503 of file G4VEmProcess.hh.

505{
506 return modelManager->SelectModel(kinEnergy, idxCouple);
507}

◆ SetBuildTableFlag()

void G4VEmProcess::SetBuildTableFlag ( G4bool val)
inline

◆ SetCrossSectionBiasingFactor()

void G4VEmProcess::SetCrossSectionBiasingFactor ( G4double f,
G4bool flag = true )

Definition at line 758 of file G4VEmProcess.cc.

759{
760 if(f > 0.0) {
761 biasFactor = f;
762 weightFlag = flag;
763 if(1 < verboseLevel) {
764 G4cout << "### SetCrossSectionBiasingFactor: for "
765 << particle->GetParticleName()
766 << " and process " << GetProcessName()
767 << " biasFactor= " << f << " weightFlag= " << flag
768 << G4endl;
769 }
770 }
771}

Referenced by G4EmExtraParameters::DefineRegParamForEM().

◆ SetCrossSectionType()

◆ SetEmMasterProcess()

void G4VEmProcess::SetEmMasterProcess ( const G4VEmProcess * ptr)
inline

Definition at line 790 of file G4VEmProcess.hh.

791{
792 masterProc = ptr;
793}

◆ SetEmModel()

void G4VEmProcess::SetEmModel ( G4VEmModel * ptr,
G4int index = 0 )

Definition at line 135 of file G4VEmProcess.cc.

136{
137 if(nullptr == ptr) { return; }
138 if(!emModels.empty()) {
139 for(auto & em : emModels) { if(em == ptr) { return; } }
140 }
141 emModels.push_back(ptr);
142}

Referenced by G4EmBuilder::ConstructElectronSSProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4EmDNABuilder::FindOrBuildAttachment(), G4EmDNABuilder::FindOrBuildChargeDecrease(), G4EmDNABuilder::FindOrBuildChargeIncrease(), G4EmDNABuilder::FindOrBuildElastic(), G4EmDNABuilder::FindOrBuildElectronSolvation(), G4EmDNABuilder::FindOrBuildExcitation(), G4EmDNABuilder::FindOrBuildIonisation(), G4EmDNABuilder::FindOrBuildVibExcitation(), G4PolarizedAnnihilation::G4PolarizedAnnihilation(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNADoubleIonisation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAElectronSolvation::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPlasmonExcitation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNAQuadrupleIonisation::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNATripleIonisation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4eeToHadrons::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4JAEAElasticScattering::InitialiseProcess(), G4MicroElecElastic::InitialiseProcess(), G4MicroElecInelastic::InitialiseProcess(), G4MicroElecLOPhononScattering::InitialiseProcess(), G4NuclearStopping::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4RayleighScattering::InitialiseProcess(), and G4EmTableUtil::PrepareEmProcess().

◆ SetEnergyOfCrossSectionMax()

void G4VEmProcess::SetEnergyOfCrossSectionMax ( std::vector< G4double > * ptr)
inline

Definition at line 684 of file G4VEmProcess.hh.

685{
687}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ SetLambdaBinning()

void G4VEmProcess::SetLambdaBinning ( G4int nbins)

Definition at line 821 of file G4VEmProcess.cc.

822{
823 if(5 < n && n < 10000000) {
824 nLambdaBins = n;
825 actBinning = true;
826 } else {
827 G4double e = (G4double)n;
828 PrintWarning("SetLambdaBinning", e);
829 }
830}

Referenced by G4GammaConversion::G4GammaConversion(), and G4PolarizedGammaConversion::G4PolarizedGammaConversion().

◆ SetLambdaTable()

void G4VEmProcess::SetLambdaTable ( G4PhysicsTable * ptr)
inline

Definition at line 662 of file G4VEmProcess.hh.

663{
664 theLambdaTable = ptr;
665}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ SetLambdaTablePrim()

void G4VEmProcess::SetLambdaTablePrim ( G4PhysicsTable * ptr)
inline

Definition at line 669 of file G4VEmProcess.hh.

670{
671 theLambdaTablePrim = ptr;
672}

Referenced by G4EmTableUtil::BuildEmProcess().

◆ SetMaxKinEnergy()

void G4VEmProcess::SetMaxKinEnergy ( G4double e)

Definition at line 846 of file G4VEmProcess.cc.

847{
848 if(minKinEnergy < e && e < 1.e+6*TeV) {
849 nLambdaBins = G4lrint(nLambdaBins*G4Log(e/minKinEnergy)
850 /G4Log(maxKinEnergy/minKinEnergy));
851 maxKinEnergy = e;
852 actMaxKinEnergy = true;
853 } else { PrintWarning("SetMaxKinEnergy", e); }
854}

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), and G4EmDNABuilder::FindOrBuildNuclearStopping().

◆ SetMinKinEnergy()

◆ SetMinKinEnergyPrim()

void G4VEmProcess::SetMinKinEnergyPrim ( G4double e)

Definition at line 858 of file G4VEmProcess.cc.

859{
860 if(theParameters->MinKinEnergy() <= e &&
861 e <= theParameters->MaxKinEnergy()) { minKinEnergyPrim = e; }
862 else { PrintWarning("SetMinKinEnergyPrim", e); }
863}
G4double MaxKinEnergy() const

Referenced by G4ComptonScattering::G4ComptonScattering(), G4JAEAElasticScattering::G4JAEAElasticScattering(), G4PhotoElectricEffect::G4PhotoElectricEffect(), G4PolarizedCompton::G4PolarizedCompton(), and G4RayleighScattering::G4RayleighScattering().

◆ SetParticle()

void G4VEmProcess::SetParticle ( const G4ParticleDefinition * p)
inlineprotected

Definition at line 733 of file G4VEmProcess.hh.

734{
735 particle = p;
736 currentParticle = p;
737}

Referenced by G4GammaGeneralProcess::G4GammaGeneralProcess(), G4eeToHadrons::InitialiseProcess(), G4GammaGeneralProcess::PreparePhysicsTable(), and PreparePhysicsTable().

◆ SetSecondaryParticle()

◆ SetSplineFlag()

◆ SetStartFromNullFlag()

◆ StartTracking()

void G4VEmProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 337 of file G4VEmProcess.cc.

338{
339 // reset parameters for the new track
340 currentParticle = track->GetParticleDefinition();
343 preStepLambda = 0.0;
344
345 if(isIon) { massRatio = proton_mass_c2/currentParticle->GetPDGMass(); }
346
347 // forced biasing only for primary particles
348 if(biasManager) {
349 if(0 == track->GetParentID()) {
350 // primary particle
351 biasFlag = true;
352 biasManager->ResetForcedInteraction();
353 }
354 }
355 for (G4int i=0; i<numberOfModels; ++i) {
356 auto ptr = GetModelByIndex(i);
357 ptr->StartTracking(track);
358 }
359}
const G4ParticleDefinition * GetParticleDefinition() const
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const

◆ StorePhysicsTable()

G4bool G4VEmProcess::StorePhysicsTable ( const G4ParticleDefinition * part,
const G4String & directory,
G4bool ascii = false )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 653 of file G4VEmProcess.cc.

656{
657 if(!isTheMaster || part != particle) { return true; }
658 if(G4EmTableUtil::StoreTable(this, part, theLambdaTable,
659 directory, "Lambda",
660 verboseLevel, ascii) &&
661 G4EmTableUtil::StoreTable(this, part, theLambdaTablePrim,
662 directory, "LambdaPrim",
663 verboseLevel, ascii)) {
664 return true;
665 }
666 return false;
667}
static G4bool StoreTable(G4VProcess *, const G4ParticleDefinition *, G4PhysicsTable *, const G4String &dir, const G4String &tname, G4int verb, G4bool ascii)

◆ StreamInfo()

void G4VEmProcess::StreamInfo ( std::ostream & outFile,
const G4ParticleDefinition & part,
G4bool rst = false ) const

Definition at line 263 of file G4VEmProcess.cc.

265{
266 G4String indent = (rst ? " " : "");
267 out << std::setprecision(6);
268 out << G4endl << indent << GetProcessName() << ": ";
269 if (!rst) {
270 out << " for " << part.GetParticleName();
271 }
272 if(fXSType != fEmNoIntegral) { out << " XStype:" << fXSType; }
273 if(applyCuts) { out << " applyCuts:1 "; }
274 G4int subtype = GetProcessSubType();
275 out << " SubType=" << subtype;
276 if (subtype == fAnnihilation) {
277 G4int mod = theParameters->PositronAtRestModelType();
278 const G4String namp[2] = {"Simple", "Allison"};
279 out << " AtRestModel:" << namp[mod];
280 }
281 if(biasFactor != 1.0) { out << " BiasingFactor=" << biasFactor; }
282 out << " BuildTable=" << buildLambdaTable << G4endl;
283 if(buildLambdaTable) {
284 if(particle == &part) {
285 for(auto & v : *theLambdaTable) {
286 if(nullptr != v) {
287 out << " Lambda table from ";
288 G4double emin = v->Energy(0);
289 G4double emax = v->GetMaxEnergy();
290 G4int nbin = G4int(v->GetVectorLength() - 1);
291 if(emin > minKinEnergy) { out << "threshold "; }
292 else { out << G4BestUnit(emin,"Energy"); }
293 out << " to "
294 << G4BestUnit(emax,"Energy")
295 << ", " << G4lrint(nbin/std::log10(emax/emin))
296 << " bins/decade, spline: "
297 << splineFlag << G4endl;
298 break;
299 }
300 }
301 } else {
302 out << " Used Lambda table of "
303 << particle->GetParticleName() << G4endl;
304 }
305 }
306 if(minKinEnergyPrim < maxKinEnergy) {
307 if(particle == &part) {
308 for(auto & v : *theLambdaTablePrim) {
309 if(nullptr != v) {
310 out << " LambdaPrime table from "
311 << G4BestUnit(v->Energy(0),"Energy")
312 << " to "
313 << G4BestUnit(v->GetMaxEnergy(),"Energy")
314 << " in " << v->GetVectorLength()-1
315 << " bins " << G4endl;
316 break;
317 }
318 }
319 } else {
320 out << " Used LambdaPrime table of "
321 << particle->GetParticleName() << G4endl;
322 }
323 }
325 modelManager->DumpModelList(out, verboseLevel);
326
327 if(verboseLevel > 2 && buildLambdaTable) {
328 out << " LambdaTable address= " << theLambdaTable << G4endl;
329 if(theLambdaTable && particle == &part) {
330 out << (*theLambdaTable) << G4endl;
331 }
332 }
333}
@ fAnnihilation
#define G4BestUnit(a, b)
virtual void StreamProcessInfo(std::ostream &) const

Referenced by G4EmTableUtil::BuildEmProcess(), and ProcessDescription().

◆ StreamProcessInfo()

virtual void G4VEmProcess::StreamProcessInfo ( std::ostream & ) const
inlineprotectedvirtual

Reimplemented in G4CoulombScattering, G4eeToHadrons, and G4eplusAnnihilation.

Definition at line 92 of file G4VEmProcess.hh.

92{};

Referenced by StreamInfo().

◆ UseBaseMaterial()

G4bool G4VEmProcess::UseBaseMaterial ( ) const
inline

Definition at line 776 of file G4VEmProcess.hh.

777{
778 return baseMat;
779}

Referenced by G4EmTableUtil::BuildEmProcess().

Member Data Documentation

◆ augerID

G4int G4VEmProcess::augerID = _AugerElectron
protected

Definition at line 402 of file G4VEmProcess.hh.

Referenced by PostStepDoIt().

◆ basedCoupleIndex

◆ baseMat

◆ biasID

G4int G4VEmProcess::biasID = _EM
protected

Definition at line 403 of file G4VEmProcess.hh.

Referenced by G4eplusAnnihilation::AtRestDoIt(), and PostStepDoIt().

◆ biasManager

◆ coupleIdxLambda

std::size_t G4VEmProcess::coupleIdxLambda = 0
protected

Definition at line 407 of file G4VEmProcess.hh.

◆ currentCouple

◆ currentCoupleIndex

◆ currentMaterial

◆ fluoID

G4int G4VEmProcess::fluoID = _Fluorescence
protected

Definition at line 401 of file G4VEmProcess.hh.

Referenced by PostStepDoIt().

◆ fParticleChange

◆ idxLambda

std::size_t G4VEmProcess::idxLambda = 0
protected

Definition at line 408 of file G4VEmProcess.hh.

◆ isTheMaster

◆ mainSecondaries

G4int G4VEmProcess::mainSecondaries = 1
protected

◆ mfpKinEnergy

G4double G4VEmProcess::mfpKinEnergy = DBL_MAX
protected

◆ preStepKinEnergy

◆ preStepLambda

◆ secID

G4int G4VEmProcess::secID = _EM
protected

◆ secParticles

std::vector<G4DynamicParticle*> G4VEmProcess::secParticles
protected

Definition at line 429 of file G4VEmProcess.hh.

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4VEmProcess(), and PostStepDoIt().

◆ theEnergyOfCrossSectionMax

std::vector<G4double>* G4VEmProcess::theEnergyOfCrossSectionMax = nullptr
protected

◆ tripletID

G4int G4VEmProcess::tripletID = _TripletElectron
protected

Definition at line 404 of file G4VEmProcess.hh.

Referenced by PostStepDoIt(), and PreparePhysicsTable().


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