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

#include <G4VEnergyLossProcess.hh>

Inheritance diagram for G4VEnergyLossProcess:

Public Member Functions

 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 ~G4VEnergyLossProcess () override
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
void ProcessDescription (std::ostream &outFile) const override
void PreparePhysicsTable (const G4ParticleDefinition &) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
void StartTracking (G4Track *) override
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) 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
G4double GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
G4double MeanFreePath (const G4Track &track)
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4VEmModelSelectModelForMaterial (G4double kinEnergy, std::size_t &idxCouple) const
void AddEmModel (G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetEmModel (G4VEmModel *, G4int index=0)
std::size_t NumberOfModels () const
G4VEmModelEmModel (std::size_t index=0) const
G4VEmModelGetModelByIndex (std::size_t idx=0, G4bool ver=false) const
void SetFluctModel (G4VEmFluctuationModel *)
G4VEmFluctuationModelFluctModel () const
void SetBaseParticle (const G4ParticleDefinition *p)
const G4ParticleDefinitionParticle () const
const G4ParticleDefinitionBaseParticle () const
const G4ParticleDefinitionSecondaryParticle () const
 G4VEnergyLossProcess (G4VEnergyLossProcess &)=delete
G4VEnergyLossProcessoperator= (const G4VEnergyLossProcess &right)=delete
void ActivateSubCutoff (const G4Region *region)
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
void ActivateForcedInteraction (G4double length, const G4String &region, G4bool flag=true)
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
void SetLossFluctuations (G4bool)
void SetSpline (G4bool val)
void SetCrossSectionType (G4CrossSectionType val)
G4CrossSectionType CrossSectionType () const
void SetIonisation (G4bool val)
G4bool IsIonisationProcess () const
void SetLinearLossLimit (G4double val)
void SetStepFunction (G4double v1, G4double v2)
void SetLowestEnergyLimit (G4double)
G4int NumberOfSubCutoffRegions () const
void SetDEDXTable (G4PhysicsTable *p, G4EmTableType tType)
void SetCSDARangeTable (G4PhysicsTable *pRange)
void SetRangeTableForLoss (G4PhysicsTable *p)
void SetInverseRangeTable (G4PhysicsTable *p)
void SetLambdaTable (G4PhysicsTable *p)
void SetTwoPeaksXS (std::vector< G4TwoPeaksXS * > *)
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
void SetDEDXBinning (G4int nbins)
void SetMinKinEnergy (G4double e)
G4double MinKinEnergy () const
void SetMaxKinEnergy (G4double e)
G4double MaxKinEnergy () const
G4double CrossSectionBiasingFactor () const
G4double GetDEDX (G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetCSDADEDX (G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDX (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
G4double GetCSDARange (G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetKineticEnergy (G4double range, const G4MaterialCutsCouple *)
G4double GetLambda (G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetLambda (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
G4bool TablesAreBuilt () const
G4PhysicsTableDEDXTable () const
G4PhysicsTableDEDXunRestrictedTable () const
G4PhysicsTableIonisationTable () const
G4PhysicsTableCSDARangeTable () const
G4PhysicsTableRangeTableForLoss () const
G4PhysicsTableInverseRangeTable () const
G4PhysicsTableLambdaTable () const
std::vector< G4TwoPeaksXS * > * TwoPeaksXS () const
std::vector< G4double > * EnergyOfCrossSectionMax () const
G4bool UseBaseMaterial () const
const G4ElementGetCurrentElement () const
void SetDynamicMassCharge (G4double massratio, G4double charge2ratio)
Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
virtual ~G4VContinuousDiscreteProcess ()
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
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 InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)=0
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
std::size_t CurrentMaterialCutsCoupleIndex () const
void SelectModel (G4double kinEnergy)
void SetParticle (const G4ParticleDefinition *p)
void SetSecondaryParticle (const G4ParticleDefinition *p)
Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
G4GPILSelection GetGPILSelection () const
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

G4ParticleChangeForLoss fParticleChange
const G4MaterialcurrentMaterial = nullptr
const G4MaterialCutsCouplecurrentCouple = nullptr
G4double preStepLambda = 0.0
G4double preStepKinEnergy = 0.0
G4double preStepScaledEnergy = 0.0
G4double mfpKinEnergy = 0.0
std::size_t currentCoupleIndex = 0
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 82 of file G4VEnergyLossProcess.hh.

Constructor & Destructor Documentation

◆ G4VEnergyLossProcess() [1/2]

G4VEnergyLossProcess::G4VEnergyLossProcess ( const G4String & name = "EnergyLoss",
G4ProcessType type = fElectromagnetic )

Definition at line 93 of file G4VEnergyLossProcess.cc.

94 :
96{
97 theParameters = G4EmParameters::Instance();
99
100 // low energy limit
101 lowestKinEnergy = theParameters->LowestElectronEnergy();
102
103 // Size of tables
104 minKinEnergy = 0.1*CLHEP::keV;
105 maxKinEnergy = 100.0*CLHEP::TeV;
106 maxKinEnergyCSDA = 1.0*CLHEP::GeV;
107 nBins = 84;
108 nBinsCSDA = 35;
109
110 invLambdaFactor = 1.0/lambdaFactor;
111
112 // default linear loss limit
113 finalRange = 1.*CLHEP::mm;
114
115 // run time objects
117 fParticleChange.SetSecondaryWeightByProcess(true);
118 modelManager = new G4EmModelManager();
120 ->GetSafetyHelper();
121 aGPILSelection = CandidateForSelection;
122
123 // initialise model
124 lManager = G4LossTableManager::Instance();
125 lManager->Register(this);
126 isMaster = lManager->IsMaster();
127
128 G4LossTableBuilder* bld = lManager->GetTableBuilder();
129 theDensityFactor = bld->GetDensityFactors();
130 theDensityIdx = bld->GetCoupleIndexes();
131 theFluctuationFlags = bld->GetFluctuationFlags();
132
133 scTracks.reserve(10);
134 secParticles.reserve(12);
135 emModels = new std::vector<G4VEmModel*>;
136}
@ CandidateForSelection
static G4EmParameters * Instance()
static const std::vector< G4double > * GetDensityFactors()
static const std::vector< G4bool > * GetFluctuationFlags()
static const std::vector< G4int > * GetCoupleIndexes()
static G4LossTableManager * Instance()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
G4ParticleChangeForLoss fParticleChange
void SetVerboseLevel(G4int value)
G4VParticleChange * pParticleChange

Referenced by BuildPhysicsTable(), G4eBremsstrahlung::G4eBremsstrahlung(), G4eIonisation::G4eIonisation(), G4ePairProduction::G4ePairProduction(), G4hhIonisation::G4hhIonisation(), G4hIonisation::G4hIonisation(), G4ionIonisation::G4ionIonisation(), G4mplIonisation::G4mplIonisation(), G4MuBremsstrahlung::G4MuBremsstrahlung(), G4MuIonisation::G4MuIonisation(), G4MuPairProduction::G4MuPairProduction(), G4PolarizedIonisation::G4PolarizedIonisation(), G4VEnergyLossProcess(), and operator=().

◆ ~G4VEnergyLossProcess()

G4VEnergyLossProcess::~G4VEnergyLossProcess ( )
override

Definition at line 140 of file G4VEnergyLossProcess.cc.

141{
142 if (isMaster) {
143 delete theEnergyOfCrossSectionMax;
144 if(nullptr != fXSpeaks) {
145 for(auto const & v : *fXSpeaks) { delete v; }
146 delete fXSpeaks;
147 }
148 }
149 delete modelManager;
150 delete biasManager;
151 delete scoffRegions;
152 delete emModels;
153 lManager->DeRegister(this);
154}

◆ G4VEnergyLossProcess() [2/2]

G4VEnergyLossProcess::G4VEnergyLossProcess ( G4VEnergyLossProcess & )
delete

Member Function Documentation

◆ ActivateForcedInteraction()

void G4VEnergyLossProcess::ActivateForcedInteraction ( G4double length,
const G4String & region,
G4bool flag = true )

Definition at line 1367 of file G4VEnergyLossProcess.cc.

1370{
1371 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
1372 if(1 < verboseLevel) {
1373 G4cout << "### ActivateForcedInteraction: for "
1374 << " process " << GetProcessName()
1375 << " length(mm)= " << length/mm
1376 << " in G4Region <" << region
1377 << "> weightFlag= " << flag
1378 << G4endl;
1379 }
1380 weightFlag = flag;
1381 biasManager->ActivateForcedInteraction(length, region);
1382}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
const G4String & GetProcessName() const

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ ActivateSecondaryBiasing()

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

Definition at line 1387 of file G4VEnergyLossProcess.cc.

1390{
1391 if (0.0 <= factor) {
1392 // Range cut can be applied only for e-
1393 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1394 { return; }
1395
1396 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
1397 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1398 if(1 < verboseLevel) {
1399 G4cout << "### ActivateSecondaryBiasing: for "
1400 << " process " << GetProcessName()
1401 << " factor= " << factor
1402 << " in G4Region <" << region
1403 << "> energyLimit(MeV)= " << energyLimit/MeV
1404 << G4endl;
1405 }
1406 }
1407}
static G4Electron * Electron()
Definition G4Electron.cc:91

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ ActivateSubCutoff()

void G4VEnergyLossProcess::ActivateSubCutoff ( const G4Region * region)

Definition at line 505 of file G4VEnergyLossProcess.cc.

506{
507 if(nullptr == scoffRegions) {
508 scoffRegions = new std::vector<const G4Region*>;
509 }
510 // the region is in the list
511 if(!scoffRegions->empty()) {
512 for (auto & reg : *scoffRegions) {
513 if (reg == r) { return; }
514 }
515 }
516 // new region
517 scoffRegions->push_back(r);
518 ++nSCoffRegions;
519}

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ AddEmModel()

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

◆ AlongStepDoIt()

G4VParticleChange * G4VEnergyLossProcess::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 763 of file G4VEnergyLossProcess.cc.

765{
766 fParticleChange.InitializeForAlongStep(track);
767 // The process has range table - calculate energy loss
768 if(!isIonisation || !currentModel->IsActive(preStepScaledEnergy)) {
769 return &fParticleChange;
770 }
771
772 G4double length = step.GetStepLength();
773 G4double eloss = 0.0;
774
775 /*
776 if(-1 < verboseLevel) {
777 const G4ParticleDefinition* d = track.GetParticleDefinition();
778 G4cout << "AlongStepDoIt for "
779 << GetProcessName() << " and particle " << d->GetParticleName()
780 << " eScaled(MeV)=" << preStepScaledEnergy/MeV
781 << " range(mm)=" << fRange/mm << " s(mm)=" << length/mm
782 << " rf=" << reduceFactor << " q^2=" << chargeSqRatio
783 << " md=" << d->GetPDGMass() << " status=" << track.GetTrackStatus()
784 << " " << track.GetMaterial()->GetName() << G4endl;
785 }
786 */
787 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
788
789 // define new weight for primary and secondaries
790 G4double weight = fParticleChange.GetParentWeight();
791 if(weightFlag) {
792 weight /= biasFactor;
793 fParticleChange.ProposeWeight(weight);
794 }
795
796 // stopping, check actual range and kinetic energy
797 if (length >= fRange || preStepKinEnergy <= lowestKinEnergy) {
798 eloss = preStepKinEnergy;
799 if (useDeexcitation) {
800 atomDeexcitation->AlongStepDeexcitation(scTracks, step,
801 eloss, (G4int)currentCoupleIndex);
802 if (!scTracks.empty()) { FillSecondariesAlongStep(weight); }
803 eloss = std::max(eloss, 0.0);
804 }
805 fParticleChange.SetProposedKineticEnergy(0.0);
806 fParticleChange.ProposeLocalEnergyDeposit(eloss);
807 return &fParticleChange;
808 }
809 // zero step length with non-zero range
810 if(length <= 0.0) { return &fParticleChange; }
811
812 // Short step
813 eloss = length*GetDEDXForScaledEnergy(preStepScaledEnergy,
814 LogScaledEkin(track));
815 /*
816 G4cout << "##### Short STEP: eloss= " << eloss
817 << " Escaled=" << preStepScaledEnergy
818 << " R=" << fRange
819 << " L=" << length
820 << " fFactor=" << fFactor << " minE=" << minKinEnergy
821 << " idxBase=" << basedCoupleIndex << G4endl;
822 */
823 // Long step
824 if(eloss > preStepKinEnergy*linLossLimit) {
825
826 const G4double x = (fRange - length)/reduceFactor;
827 const G4double de = preStepKinEnergy - ScaledKinEnergyForLoss(x)/massRatio;
828 if(de > 0.0) { eloss = de; }
829 /*
830 if(-1 < verboseLevel)
831 G4cout << " Long STEP: rPre(mm)="
832 << GetScaledRangeForScaledEnergy(preStepScaledEnergy)/mm
833 << " x(mm)=" << x/mm
834 << " eloss(MeV)=" << eloss/MeV
835 << " rFactor=" << reduceFactor
836 << " massRatio=" << massRatio
837 << G4endl;
838 */
839 }
840
841 /*
842 if(-1 < verboseLevel ) {
843 G4cout << "Before fluct: eloss(MeV)= " << eloss/MeV
844 << " e-eloss= " << preStepKinEnergy-eloss
845 << " step(mm)= " << length/mm << " range(mm)= " << fRange/mm
846 << " fluct= " << lossFluctuationFlag << G4endl;
847 }
848 */
849
850 const G4double cut = (*theCuts)[currentCoupleIndex];
851 G4double esec = 0.0;
852
853 // Corrections, which cannot be tabulated
854 if (isIon) {
855 currentModel->SetCurrentCouple(currentCouple);
856 currentModel->CorrectionsAlongStep(currentMaterial,
857 track.GetParticleDefinition(),
858 preStepKinEnergy, cut,
859 length, eloss);
860 eloss = std::max(eloss, 0.0);
861 }
862
863 // Sample fluctuations if not full energy loss
864 if(eloss >= preStepKinEnergy) {
865 eloss = preStepKinEnergy;
866
867 } else if ((*theFluctuationFlags)[currentCoupleIndex]) {
868 const G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
869 const G4double tcut = std::min(cut, tmax);
870 G4VEmFluctuationModel* fluc = currentModel->GetModelOfFluctuations();
871 eloss = fluc->SampleFluctuations(currentCouple,dynParticle,
872 tcut, tmax, length, eloss);
873 /*
874 if(-1 < verboseLevel)
875 G4cout << "After fluct: eloss(MeV)= " << eloss/MeV
876 << " fluc= " << (eloss-eloss0)/MeV
877 << " ChargeSqRatio= " << chargeSqRatio
878 << " massRatio= " << massRatio << " tmax= " << tmax << G4endl;
879 */
880 }
881
882 // deexcitation
883 if (useDeexcitation) {
884 G4double esecfluo = preStepKinEnergy;
885 G4double de = esecfluo;
886 atomDeexcitation->AlongStepDeexcitation(scTracks, step,
888
889 // sum of de-excitation energies
890 esecfluo -= de;
891
892 // subtracted from energy loss
893 if(eloss >= esecfluo) {
894 esec += esecfluo;
895 eloss -= esecfluo;
896 } else {
897 esec += esecfluo;
898 eloss = 0.0;
899 }
900 }
901 if(nullptr != subcutProducer && IsRegionForCubcutProcessor(track)) {
902 subcutProducer->SampleSecondaries(step, scTracks, eloss, cut);
903 }
904 // secondaries from atomic de-excitation and subcut
905 if(!scTracks.empty()) { FillSecondariesAlongStep(weight); }
906
907 // Energy balance
908 G4double finalT = preStepKinEnergy - eloss - esec;
909 if (finalT <= lowestKinEnergy) {
910 eloss += finalT;
911 finalT = 0.0;
912 } else if(isIon) {
913 fParticleChange.SetProposedCharge(
914 currentModel->GetParticleCharge(track.GetParticleDefinition(),
915 currentMaterial,finalT));
916 }
917 eloss = std::max(eloss, 0.0);
918
919 fParticleChange.SetProposedKineticEnergy(finalT);
920 fParticleChange.ProposeLocalEnergyDeposit(eloss);
921 /*
922 if(-1 < verboseLevel) {
923 G4double del = finalT + eloss + esec - preStepKinEnergy;
924 G4cout << "Final value eloss(MeV)= " << eloss/MeV
925 << " preStepKinEnergy= " << preStepKinEnergy
926 << " postStepKinEnergy= " << finalT
927 << " de(keV)= " << del/keV
928 << " lossFlag= " << lossFluctuationFlag
929 << " status= " << track.GetTrackStatus()
930 << G4endl;
931 }
932 */
933 return &fParticleChange;
934}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetStepLength() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VEnergyLossProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 561 of file G4VEnergyLossProcess.cc.

564{
565 G4double x = DBL_MAX;
566 *selection = aGPILSelection;
567 if(isIonisation && currentModel->IsActive(preStepScaledEnergy)) {
568 GetScaledRangeForScaledEnergy(preStepScaledEnergy, LogScaledEkin(track));
569 x = (useCutAsFinalRange) ? std::min(finalRange,
570 currentCouple->GetProductionCuts()->GetProductionCut(1)) : finalRange;
571 x = (fRange > x) ? fRange*dRoverRange + x*(1.0 - dRoverRange)*(2.0 - x/fRange)
572 : fRange;
573 /*
574 G4cout<<"AlongStepGPIL: " << GetProcessName()<<": e="<<preStepKinEnergy
575 << " fRange=" << fRange << " finR=" << finR <<" stepLimit="<<x<<G4endl;
576 */
577 }
578 return x;
579}
#define DBL_MAX
Definition templates.hh:62

Referenced by ContinuousStepLimit().

◆ BaseParticle()

const G4ParticleDefinition * G4VEnergyLossProcess::BaseParticle ( ) const
inline

◆ BuildDEDXTable()

G4PhysicsTable * G4VEnergyLossProcess::BuildDEDXTable ( G4EmTableType tType = fRestricted)

Definition at line 405 of file G4VEnergyLossProcess.cc.

406{
407 G4PhysicsTable* table = nullptr;
408 G4double emax = maxKinEnergy;
409 G4int bin = nBins;
410
411 if(fTotal == tType) {
412 emax = maxKinEnergyCSDA;
413 bin = nBinsCSDA;
414 table = theDEDXunRestrictedTable;
415 } else if(fRestricted == tType) {
416 table = theDEDXTable;
417 } else {
418 G4cout << "G4VEnergyLossProcess::BuildDEDXTable WARNING: wrong type "
419 << tType << G4endl;
420 }
421 if(1 < verboseLevel) {
422 G4cout << "G4VEnergyLossProcess::BuildDEDXTable() of type " << tType
423 << " for " << GetProcessName()
424 << " and " << particle->GetParticleName()
425 << "spline=" << spline << G4endl;
426 }
427 if(nullptr == table) { return table; }
428
429 G4LossTableBuilder* bld = lManager->GetTableBuilder();
430 G4EmTableUtil::BuildDEDXTable(this, particle, modelManager, bld,
431 table, minKinEnergy, emax, bin,
432 verboseLevel, tType, spline);
433 return table;
434}
@ fTotal
@ fRestricted
static void BuildDEDXTable(G4VEnergyLossProcess *proc, const G4ParticleDefinition *part, G4EmModelManager *modelManager, G4LossTableBuilder *bld, G4PhysicsTable *table, const G4double minKinEnergy, const G4double maxKinEnergy, const G4int nbins, const G4int verbose, const G4EmTableType tType, const G4bool splineFlag)

◆ BuildLambdaTable()

G4PhysicsTable * G4VEnergyLossProcess::BuildLambdaTable ( G4EmTableType tType = fRestricted)

Definition at line 438 of file G4VEnergyLossProcess.cc.

439{
440 if(nullptr == theLambdaTable) { return theLambdaTable; }
441
442 G4double scale = theParameters->MaxKinEnergy()/theParameters->MinKinEnergy();
443 G4int nbin =
444 theParameters->NumberOfBinsPerDecade()*G4lrint(std::log10(scale));
445 scale = nbin/G4Log(scale);
446
447 G4LossTableBuilder* bld = lManager->GetTableBuilder();
448 G4EmTableUtil::BuildLambdaTable(this, particle, modelManager,
449 bld, theLambdaTable, theCuts,
450 minKinEnergy, maxKinEnergy, scale,
451 verboseLevel, spline);
452 return theLambdaTable;
453}
G4double G4Log(G4double x)
Definition G4Log.hh:169
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

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 336 of file G4VEnergyLossProcess.cc.

337{
338 if(1 < verboseLevel) {
339 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() for "
340 << GetProcessName()
341 << " and particle " << part.GetParticleName()
342 << "; the first particle " << particle->GetParticleName();
343 if(baseParticle) {
344 G4cout << "; base: " << baseParticle->GetParticleName();
345 }
346 G4cout << G4endl;
347 G4cout << " TablesAreBuilt= " << tablesAreBuilt << " isIon= " << isIon
348 << " spline=" << spline << " ptr: " << this << G4endl;
349 }
350
351 if(&part == particle) {
352 if(isMaster) {
353 lManager->BuildPhysicsTable(particle, this);
354
355 } else {
356 const auto masterProcess =
357 static_cast<const G4VEnergyLossProcess*>(GetMasterProcess());
358
359 numberOfModels = modelManager->NumberOfModels();
360 G4EmTableUtil::BuildLocalElossProcess(this, masterProcess,
361 particle, numberOfModels);
362 tablesAreBuilt = true;
363 baseMat = masterProcess->UseBaseMaterial();
364 lManager->LocalPhysicsTables(particle, this);
365 }
366
367 // needs to be done only once
368 safetyHelper->InitialiseHelper();
369 }
370 // Added tracking cut to avoid tracking artifacts
371 // and identified deexcitation flag
372 if(isIonisation) {
373 atomDeexcitation = lManager->AtomDeexcitation();
374 if(nullptr != atomDeexcitation) {
375 if(atomDeexcitation->IsPIXEActive()) { useDeexcitation = true; }
376 }
377 }
378
379 // protection against double printout
380 if(theParameters->IsPrintLocked()) { return; }
381
382 // explicitly defined printout by particle name
383 G4String num = part.GetParticleName();
384 if(1 < verboseLevel ||
385 (0 < verboseLevel && (num == "e-" ||
386 num == "e+" || num == "mu+" ||
387 num == "mu-" || num == "proton"||
388 num == "pi+" || num == "pi-" ||
389 num == "kaon+" || num == "kaon-" ||
390 num == "alpha" || num == "anti_proton" ||
391 num == "GenericIon"|| num == "alpha+" ))) {
392 StreamInfo(G4cout, part);
393 }
394 if(1 < verboseLevel) {
395 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() done for "
396 << GetProcessName()
397 << " and particle " << part.GetParticleName();
398 if(isIonisation) { G4cout << " isIonisation flag=1"; }
399 G4cout << " baseMat=" << baseMat << G4endl;
400 }
401}
static void BuildLocalElossProcess(G4VEnergyLossProcess *proc, const G4VEnergyLossProcess *masterProc, const G4ParticleDefinition *part, const G4int nModels)
const G4String & GetParticleName() const
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
const G4VProcess * GetMasterProcess() const

Referenced by G4PolarizedIonisation::BuildPhysicsTable().

◆ ContinuousStepLimit()

G4double G4VEnergyLossProcess::ContinuousStepLimit ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )

Definition at line 1208 of file G4VEnergyLossProcess.cc.

1211{
1212 return AlongStepGetPhysicalInteractionLength(track, x, y, z, &aGPILSelection);
1213}
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override

◆ CrossSectionBiasingFactor()

G4double G4VEnergyLossProcess::CrossSectionBiasingFactor ( ) const
inline

Definition at line 916 of file G4VEnergyLossProcess.hh.

917{
918 return biasFactor;
919}

◆ CrossSectionPerVolume() [1/2]

G4double G4VEnergyLossProcess::CrossSectionPerVolume ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple )

◆ CrossSectionPerVolume() [2/2]

G4double G4VEnergyLossProcess::CrossSectionPerVolume ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )

Definition at line 1175 of file G4VEnergyLossProcess.cc.

1178{
1179 // Cross section per volume is calculated
1180 DefineMaterial(couple);
1181 G4double cross = 0.0;
1182 if (nullptr != theLambdaTable) {
1183 cross = GetLambdaForScaledEnergy(kineticEnergy * massRatio,
1184 logKineticEnergy + logMassRatio);
1185 } else {
1186 SelectModel(kineticEnergy*massRatio);
1187 cross = (!baseMat) ? biasFactor : biasFactor*(*theDensityFactor)[currentCoupleIndex];
1188 cross *= (currentModel->CrossSectionPerVolume(currentMaterial, particle, kineticEnergy,
1189 (*theCuts)[currentCoupleIndex]));
1190 }
1191 return std::max(cross, 0.0);
1192}
void SelectModel(G4double kinEnergy)

◆ CrossSectionType()

G4CrossSectionType G4VEnergyLossProcess::CrossSectionType ( ) const
inline

Definition at line 881 of file G4VEnergyLossProcess.hh.

882{
883 return fXSType;
884}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ CSDARangeTable()

G4PhysicsTable * G4VEnergyLossProcess::CSDARangeTable ( ) const
inline

Definition at line 951 of file G4VEnergyLossProcess.hh.

952{
953 return theCSDARangeTable;
954}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ CurrentMaterialCutsCoupleIndex()

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

Definition at line 537 of file G4VEnergyLossProcess.hh.

538{
539 return currentCoupleIndex;
540}

◆ DEDXTable()

G4PhysicsTable * G4VEnergyLossProcess::DEDXTable ( ) const
inline

◆ DEDXunRestrictedTable()

G4PhysicsTable * G4VEnergyLossProcess::DEDXunRestrictedTable ( ) const
inline

Definition at line 937 of file G4VEnergyLossProcess.hh.

938{
939 return theDEDXunRestrictedTable;
940}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ EmModel()

◆ EnergyOfCrossSectionMax()

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

Definition at line 987 of file G4VEnergyLossProcess.hh.

988{
989 return theEnergyOfCrossSectionMax;
990}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ FluctModel()

◆ GetContinuousStepLimit()

G4double G4VEnergyLossProcess::GetContinuousStepLimit ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 1229 of file G4VEnergyLossProcess.cc.

1232{
1233 return DBL_MAX;
1234}

◆ GetCSDADEDX()

G4double G4VEnergyLossProcess::GetCSDADEDX ( G4double kineticEnergy,
const G4MaterialCutsCouple *  )
inline

◆ GetCSDARange()

G4double G4VEnergyLossProcess::GetCSDARange ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple )
inline

Definition at line 764 of file G4VEnergyLossProcess.hh.

766{
767 DefineMaterial(couple);
768 return (nullptr == theCSDARangeTable) ? DBL_MAX :
769 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
770}

◆ GetCurrentElement()

const G4Element * G4VEnergyLossProcess::GetCurrentElement ( ) const

Definition at line 1342 of file G4VEnergyLossProcess.cc.

1343{
1344 return (nullptr != currentModel)
1345 ? currentModel->GetCurrentElement(currentMaterial) : nullptr;
1346}

◆ GetDEDX() [1/2]

G4double G4VEnergyLossProcess::GetDEDX ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple )
inline

Definition at line 722 of file G4VEnergyLossProcess.hh.

724{
725 DefineMaterial(couple);
726 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
727}

◆ GetDEDX() [2/2]

G4double G4VEnergyLossProcess::GetDEDX ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )
inline

Definition at line 732 of file G4VEnergyLossProcess.hh.

735{
736 DefineMaterial(couple);
737 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
738}

◆ GetDEDXDispersion()

G4double G4VEnergyLossProcess::GetDEDXDispersion ( const G4MaterialCutsCouple * couple,
const G4DynamicParticle * dp,
G4double length )

Definition at line 1156 of file G4VEnergyLossProcess.cc.

1160{
1161 DefineMaterial(couple);
1162 G4double ekin = dp->GetKineticEnergy();
1163 SelectModel(ekin*massRatio);
1164 G4double tmax = currentModel->MaxSecondaryKinEnergy(dp);
1165 G4double tcut = std::min(tmax,(*theCuts)[currentCoupleIndex]);
1166 G4double d = 0.0;
1167 G4VEmFluctuationModel* fm = currentModel->GetModelOfFluctuations();
1168 if(nullptr != fm) { d = fm->Dispersion(currentMaterial,dp,tcut,tmax,length); }
1169 return d;
1170}
G4double GetKineticEnergy() const
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length)=0

◆ GetKineticEnergy()

G4double G4VEnergyLossProcess::GetKineticEnergy ( G4double range,
const G4MaterialCutsCouple * couple )
inline

Definition at line 775 of file G4VEnergyLossProcess.hh.

777{
778 DefineMaterial(couple);
779 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
780}

◆ GetLambda() [1/2]

G4double G4VEnergyLossProcess::GetLambda ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple )
inline

Definition at line 785 of file G4VEnergyLossProcess.hh.

787{
788 DefineMaterial(couple);
789 return (nullptr != theLambdaTable) ?
790 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
791}

◆ GetLambda() [2/2]

G4double G4VEnergyLossProcess::GetLambda ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )
inline

Definition at line 796 of file G4VEnergyLossProcess.hh.

799{
800 DefineMaterial(couple);
801 return (nullptr != theLambdaTable) ?
802 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
803 : 0.0;
804}

◆ GetMeanFreePath()

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

Implements G4VContinuousDiscreteProcess.

Definition at line 1217 of file G4VEnergyLossProcess.cc.

1222{
1224 return MeanFreePath(track);
1225}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double MeanFreePath(const G4Track &track)

Referenced by G4PolarizedIonisation::GetMeanFreePath().

◆ GetModelByIndex()

G4VEmModel * G4VEnergyLossProcess::GetModelByIndex ( std::size_t idx = 0,
G4bool ver = false ) const
inline

Definition at line 1016 of file G4VEnergyLossProcess.hh.

1017{
1018 return modelManager->GetModel((G4int)idx, ver);
1019}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ GetRange() [1/2]

G4double G4VEnergyLossProcess::GetRange ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple )
inline

Definition at line 743 of file G4VEnergyLossProcess.hh.

745{
746 DefineMaterial(couple);
747 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
748}

◆ GetRange() [2/2]

G4double G4VEnergyLossProcess::GetRange ( G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )
inline

Definition at line 753 of file G4VEnergyLossProcess.hh.

756{
757 DefineMaterial(couple);
758 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
759}

◆ InitialiseEnergyLossProcess()

◆ InverseRangeTable()

G4PhysicsTable * G4VEnergyLossProcess::InverseRangeTable ( ) const
inline

Definition at line 965 of file G4VEnergyLossProcess.hh.

966{
967 return theInverseRangeTable;
968}

Referenced by G4EmTableUtil::BuildLocalElossProcess(), G4LossTableManager::LocalPhysicsTables(), and G4EmCalculator::PrintInverseRangeTable().

◆ IonisationTable()

G4PhysicsTable * G4VEnergyLossProcess::IonisationTable ( ) const
inline

Definition at line 944 of file G4VEnergyLossProcess.hh.

945{
946 return theIonisationTable;
947}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ IsIonisationProcess()

G4bool G4VEnergyLossProcess::IsIonisationProcess ( ) const
inline

◆ LambdaPhysicsVector()

G4PhysicsVector * G4VEnergyLossProcess::LambdaPhysicsVector ( const G4MaterialCutsCouple * couple,
G4double cut )
protected

Definition at line 1239 of file G4VEnergyLossProcess.cc.

1241{
1242 DefineMaterial(couple);
1243 G4PhysicsVector* v = (*theLambdaTable)[basedCoupleIndex];
1244 return new G4PhysicsVector(*v);
1245}

◆ LambdaTable()

G4PhysicsTable * G4VEnergyLossProcess::LambdaTable ( ) const
inline

Definition at line 972 of file G4VEnergyLossProcess.hh.

973{
974 return theLambdaTable;
975}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ MaxKinEnergy()

G4double G4VEnergyLossProcess::MaxKinEnergy ( ) const
inline

Definition at line 909 of file G4VEnergyLossProcess.hh.

910{
911 return maxKinEnergy;
912}

◆ MeanFreePath()

G4double G4VEnergyLossProcess::MeanFreePath ( const G4Track & track)

Definition at line 1196 of file G4VEnergyLossProcess.cc.

1197{
1198 DefineMaterial(track.GetMaterialCutsCouple());
1199 const G4double kinEnergy = track.GetKineticEnergy();
1200 const G4double logKinEnergy = track.GetDynamicParticle()->GetLogKineticEnergy();
1201 const G4double cs = GetLambdaForScaledEnergy(kinEnergy * massRatio,
1202 logKinEnergy + logMassRatio);
1203 return (0.0 < cs) ? 1.0/cs : DBL_MAX;
1204}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

Referenced by GetMeanFreePath().

◆ MinKinEnergy()

G4double G4VEnergyLossProcess::MinKinEnergy ( ) const
inline

Definition at line 902 of file G4VEnergyLossProcess.hh.

903{
904 return minKinEnergy;
905}

◆ MinPrimaryEnergy()

G4double G4VEnergyLossProcess::MinPrimaryEnergy ( const G4ParticleDefinition * ,
const G4Material * ,
G4double cut )
virtual

◆ NumberOfModels()

std::size_t G4VEnergyLossProcess::NumberOfModels ( ) const
inline

Definition at line 1001 of file G4VEnergyLossProcess.hh.

1002{
1003 return numberOfModels;
1004}

◆ NumberOfSubCutoffRegions()

G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions ( ) const
inline

Definition at line 895 of file G4VEnergyLossProcess.hh.

896{
897 return nSCoffRegions;
898}

◆ operator=()

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

◆ Particle()

const G4ParticleDefinition * G4VEnergyLossProcess::Particle ( ) const
inline

Definition at line 845 of file G4VEnergyLossProcess.hh.

846{
847 return particle;
848}

Referenced by G4LossTableManager::BuildPhysicsTable(), and G4LossTableManager::LocalPhysicsTables().

◆ PostStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 978 of file G4VEnergyLossProcess.cc.

980{
981 // clear number of interaction lengths in any case
984
985 fParticleChange.InitializeForPostStep(track);
986 const G4double finalT = track.GetKineticEnergy();
987
988 const G4double postStepScaledEnergy = finalT*massRatio;
989 SelectModel(postStepScaledEnergy);
990
991 if(!currentModel->IsActive(postStepScaledEnergy)) {
992 return &fParticleChange;
993 }
994 /*
995 if(1 < verboseLevel) {
996 G4cout<<GetProcessName()<<" PostStepDoIt: E(MeV)= "<< finalT/MeV<< G4endl;
997 }
998 */
999 // forced process - should happen only once per track
1000 if(biasFlag) {
1001 if(biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
1002 biasFlag = false;
1003 }
1004 }
1005 const G4DynamicParticle* dp = track.GetDynamicParticle();
1006
1007 // Integral approach
1008 if (fXSType != fEmNoIntegral) {
1009 const G4double logFinalT = dp->GetLogKineticEnergy();
1010 G4double lx = GetLambdaForScaledEnergy(postStepScaledEnergy,
1011 logFinalT + logMassRatio);
1012 lx = std::max(lx, 0.0);
1013
1014 // if both lg and lx are zero then no interaction
1015 if(preStepLambda*G4UniformRand() >= lx) {
1016 return &fParticleChange;
1017 }
1018 }
1019
1020 // define new weight for primary and secondaries
1021 G4double weight = fParticleChange.GetParentWeight();
1022 if(weightFlag) {
1023 weight /= biasFactor;
1024 fParticleChange.ProposeWeight(weight);
1025 }
1026
1027 const G4double tcut = (*theCuts)[currentCoupleIndex];
1028
1029 // sample secondaries
1030 secParticles.clear();
1031 currentModel->SampleSecondaries(&secParticles, currentCouple, dp, tcut);
1032
1033 const G4int num0 = (G4int)secParticles.size();
1034
1035 // bremsstrahlung splitting or Russian roulette
1036 if(biasManager) {
1037 if(biasManager->SecondaryBiasingRegion((G4int)currentCoupleIndex)) {
1038 G4double eloss = 0.0;
1039 weight *= biasManager->ApplySecondaryBiasing(
1040 secParticles,
1041 track, currentModel,
1042 &fParticleChange, eloss,
1043 (G4int)currentCoupleIndex, tcut,
1044 step.GetPostStepPoint()->GetSafety());
1045 if(eloss > 0.0) {
1046 eloss += fParticleChange.GetLocalEnergyDeposit();
1047 fParticleChange.ProposeLocalEnergyDeposit(eloss);
1048 }
1049 }
1050 }
1051
1052 // save secondaries
1053 const G4int num = (G4int)secParticles.size();
1054 if(num > 0) {
1055
1056 fParticleChange.SetNumberOfSecondaries(num);
1057 G4double time = track.GetGlobalTime();
1058
1059 G4int n1(0), n2(0);
1060 if(num0 > mainSecondaries) {
1061 currentModel->FillNumberOfSecondaries(n1, n2);
1062 }
1063
1064 for (G4int i=0; i<num; ++i) {
1065 if(nullptr != secParticles[i]) {
1066 G4Track* t = new G4Track(secParticles[i], time, track.GetPosition());
1068 if (biasManager) {
1069 t->SetWeight(weight * biasManager->GetWeight(i));
1070 } else {
1071 t->SetWeight(weight);
1072 }
1073 if(i < num0) {
1074 t->SetCreatorModelID(secID);
1075 } else if(i < num0 + n1) {
1076 t->SetCreatorModelID(tripletID);
1077 } else {
1078 t->SetCreatorModelID(biasID);
1079 }
1080
1081 //G4cout << "Secondary(post step) has weight " << t->GetWeight()
1082 // << ", kenergy " << t->GetKineticEnergy()/MeV << " MeV"
1083 // << " time= " << time/ns << " ns " << G4endl;
1084 pParticleChange->AddSecondary(t);
1085 }
1086 }
1087 }
1088
1089 if(0.0 == fParticleChange.GetProposedKineticEnergy() &&
1090 fAlive == fParticleChange.GetTrackStatus()) {
1091 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
1092 { fParticleChange.ProposeTrackStatus(fStopButAlive); }
1093 else { fParticleChange.ProposeTrackStatus(fStopAndKill); }
1094 }
1095
1096 /*
1097 if(-1 < verboseLevel) {
1098 G4cout << "::PostStepDoIt: Sample secondary; Efin= "
1099 << fParticleChange.GetProposedKineticEnergy()/MeV
1100 << " MeV; model= (" << currentModel->LowEnergyLimit()
1101 << ", " << currentModel->HighEnergyLimit() << ")"
1102 << " preStepLambda= " << preStepLambda
1103 << " dir= " << track.GetMomentumDirection()
1104 << " status= " << track.GetTrackStatus()
1105 << G4endl;
1106 }
1107 */
1108 return &fParticleChange;
1109}
@ fEmNoIntegral
@ fAlive
@ fStopAndKill
@ fStopButAlive
#define G4UniformRand()
Definition Randomize.hh:52
G4double GetSafety() const
G4StepPoint * GetPostStepPoint() 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)
G4double theNumberOfInteractionLengthLeft

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 583 of file G4VEnergyLossProcess.cc.

587{
588 // condition is set to "Not Forced"
590 G4double x = DBL_MAX;
591
592 // initialisation of material, mass, charge, model
593 // at the beginning of the step
594 DefineMaterial(track.GetMaterialCutsCouple());
598
599 if(!currentModel->IsActive(preStepScaledEnergy)) {
602 preStepLambda = 0.0;
604 return x;
605 }
606
607 // change effective charge of a charged particle on fly
608 if(isIon) {
609 const G4double q2 = currentModel->ChargeSquareRatio(track);
610 fFactor = q2*biasFactor;
611 if(baseMat) { fFactor *= (*theDensityFactor)[currentCoupleIndex]; }
612 reduceFactor = 1.0/(fFactor*massRatio);
613 if (lossFluctuationFlag) {
614 auto fluc = currentModel->GetModelOfFluctuations();
615 fluc->SetParticleAndCharge(track.GetParticleDefinition(), q2);
616 }
617 }
618
619 // forced biasing only for primary particles
620 if(biasManager) {
621 if(0 == track.GetParentID() && biasFlag &&
622 biasManager->ForcedInteractionRegion((G4int)currentCoupleIndex)) {
623 return biasManager->GetStepLimit((G4int)currentCoupleIndex, previousStepSize);
624 }
625 }
626
627 ComputeLambdaForScaledEnergy(preStepScaledEnergy, track);
628
629 // zero cross section
630 if(preStepLambda <= 0.0) {
633 } else {
634
635 // non-zero cross section
637
638 // beggining of tracking (or just after DoIt of this process)
641
642 } else if(currentInteractionLength < DBL_MAX) {
643
644 // subtract NumberOfInteractionLengthLeft using previous step
646 previousStepSize/currentInteractionLength;
647
650 }
651
652 // new mean free path and step limit
655 }
656#ifdef G4VERBOSE
657 if (verboseLevel>2) {
658 G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ";
659 G4cout << "[ " << GetProcessName() << "]" << G4endl;
660 G4cout << " for " << track.GetParticleDefinition()->GetParticleName()
661 << " in Material " << currentMaterial->GetName()
662 << " Ekin(MeV)= " << preStepKinEnergy/MeV
663 << " track material: " << track.GetMaterial()->GetName()
664 <<G4endl;
665 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]"
666 << "InteractionLength= " << x/cm <<"[cm] " <<G4endl;
667 }
668#endif
669 return x;
670}
const G4String & GetName() const
G4Material * GetMaterial() const
G4int GetParentID() const
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength

Referenced by G4PolarizedIonisation::PostStepGetPhysicalInteractionLength().

◆ PreparePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 204 of file G4VEnergyLossProcess.cc.

205{
206 particle = G4EmTableUtil::CheckIon(this, &part, particle,
207 verboseLevel, isIon);
208
209 if( particle != &part ) {
210 if(!isIon) { lManager->RegisterExtraParticle(&part, this); }
211 if(1 < verboseLevel) {
212 G4cout << "### G4VEnergyLossProcess::PreparePhysicsTable()"
213 << " interrupted for " << GetProcessName() << " and "
214 << part.GetParticleName() << " isIon=" << isIon
215 << " spline=" << spline << G4endl;
216 }
217 return;
218 }
219
220 tablesAreBuilt = false;
221 if (GetProcessSubType() == fIonisation) { SetIonisation(true); }
222
223 G4LossTableBuilder* bld = lManager->GetTableBuilder();
224 lManager->PreparePhysicsTable(&part, this);
225
226 // Base particle and set of models can be defined here
227 InitialiseEnergyLossProcess(particle, baseParticle);
228
229 // parameters of the process
230 lossFluctuationFlag = theParameters->LossFluctuation();
231 useCutAsFinalRange = theParameters->UseCutAsFinalRange();
232 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
233 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
234 if(!actBinning) { nBins = theParameters->NumberOfBins(); }
235 maxKinEnergyCSDA = theParameters->MaxEnergyForCSDARange();
236 nBinsCSDA = theParameters->NumberOfBinsPerDecade()
237 *G4lrint(std::log10(maxKinEnergyCSDA/minKinEnergy));
238 if(!actLinLossLimit) { linLossLimit = theParameters->LinearLossLimit(); }
239 lambdaFactor = theParameters->LambdaFactor();
240 invLambdaFactor = 1.0/lambdaFactor;
241 if(isMaster) { SetVerboseLevel(theParameters->Verbose()); }
242 else { SetVerboseLevel(theParameters->WorkerVerbose()); }
243 // integral option may be disabled
244 if(!theParameters->Integral()) { fXSType = fEmNoIntegral; }
245
246 theParameters->DefineRegParamForLoss(this);
247
248 fRangeEnergy = 0.0;
249
250 G4double initialCharge = particle->GetPDGCharge();
251 G4double initialMass = particle->GetPDGMass();
252
253 theParameters->FillStepFunction(particle, this);
254
255 // parameters for scaling from the base particle
256 if (nullptr != baseParticle) {
257 massRatio = (baseParticle->GetPDGMass())/initialMass;
258 logMassRatio = G4Log(massRatio);
259 G4double q = initialCharge/baseParticle->GetPDGCharge();
260 chargeSqRatio = q*q;
261 if(chargeSqRatio > 0.0) { reduceFactor = 1.0/(chargeSqRatio*massRatio); }
262 }
263 lowestKinEnergy = (initialMass < CLHEP::MeV)
264 ? theParameters->LowestElectronEnergy()
265 : theParameters->LowestMuHadEnergy();
266
267 // Tables preparation
268 if (isMaster && nullptr == baseParticle) {
269 if (nullptr == theData) { theData = new G4EmDataHandler(7, particle->GetParticleName()); }
270
271 if (isIonisation && theDEDXTable != theIonisationTable) {
272 theData->CleanTable(0);
273 theDEDXTable = theIonisationTable;
274 theIonisationTable = nullptr;
275 }
276 theDEDXTable = theData->MakeTable(theDEDXTable, 0);
277 bld->InitialiseBaseMaterials(theDEDXTable);
278 theData->MakeTable(theIonisationTable, 1);
279
280 if (theParameters->BuildCSDARange()) {
281 theDEDXunRestrictedTable = theData->MakeTable(theDEDXunRestrictedTable, 2);
282 if(isIonisation) { theCSDARangeTable = theData->MakeTable(theCSDARangeTable, 3); }
283 }
284
285 theLambdaTable = theData->MakeTable(theLambdaTable, 4);
286 if(isIonisation) {
287 theRangeTableForLoss = theData->MakeTable(theRangeTableForLoss, 5);
288 theInverseRangeTable = theData->MakeTable(theInverseRangeTable, 6);
289 }
290 }
291
292 // forced biasing
293 if(nullptr != biasManager) {
294 biasManager->Initialise(part,GetProcessName(),verboseLevel);
295 biasFlag = false;
296 }
297 baseMat = bld->GetBaseMaterialFlag();
298 numberOfModels = modelManager->NumberOfModels();
299 currentModel = modelManager->GetModel(0);
300 G4EmTableUtil::UpdateModels(this, modelManager, maxKinEnergy,
301 numberOfModels, secID, biasID,
302 mainSecondaries, baseMat, isMaster,
303 theParameters->UseAngularGeneratorForIonisation());
304 theCuts = modelManager->Initialise(particle, secondaryParticle,
306 // subcut processor
307 if(isIonisation) {
308 subcutProducer = lManager->SubCutProducer();
309 }
310
311 if(1 < verboseLevel) {
312 G4cout << "G4VEnergyLossProcess::PrepearPhysicsTable() is done "
313 << " for " << GetProcessName() << " and " << particle->GetParticleName()
314 << " isIon= " << isIon << " spline=" << spline;
315 if(baseParticle) {
316 G4cout << "; base: " << baseParticle->GetParticleName();
317 }
318 G4cout << G4endl;
319 G4cout << " chargeSqRatio= " << chargeSqRatio
320 << " massRatio= " << massRatio
321 << " reduceFactor= " << reduceFactor << G4endl;
322 if (nSCoffRegions > 0) {
323 G4cout << " SubCut secondary production is ON for regions: " << G4endl;
324 for (G4int i=0; i<nSCoffRegions; ++i) {
325 const G4Region* r = (*scoffRegions)[i];
326 G4cout << " " << r->GetName() << G4endl;
327 }
328 } else if(nullptr != subcutProducer) {
329 G4cout << " SubCut secondary production is ON for all regions" << G4endl;
330 }
331 }
332}
@ fIonisation
static void UpdateModels(G4VEnergyLossProcess *proc, G4EmModelManager *modelManager, const G4double maxKinEnergy, const G4int nModels, G4int &secID, G4int &biasID, G4int &mainSecondaries, const G4bool baseMat, const G4bool isMaster, const G4bool useAGen)
static const G4ParticleDefinition * CheckIon(G4VEnergyLossProcess *proc, const G4ParticleDefinition *part, const G4ParticleDefinition *particle, const G4int verboseLevel, G4bool &isIon)
static G4bool GetBaseMaterialFlag()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
const G4String & GetName() const
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
G4int GetProcessSubType() const

◆ ProcessDescription()

◆ RangeTableForLoss()

G4PhysicsTable * G4VEnergyLossProcess::RangeTableForLoss ( ) const
inline

Definition at line 958 of file G4VEnergyLossProcess.hh.

959{
960 return theRangeTableForLoss;
961}

Referenced by G4EmTableUtil::BuildLocalElossProcess(), G4LossTableManager::LocalPhysicsTables(), and G4EmCalculator::PrintRangeTable().

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 1139 of file G4VEnergyLossProcess.cc.

1141{
1142 if (!isMaster || nullptr != baseParticle || part != particle ) return true;
1143 for(std::size_t i=0; i<7; ++i) {
1144 // ionisation table only for ionisation process
1145 if (!isIonisation && 1 == i) { continue; }
1146 if(!G4EmTableUtil::RetrieveTable(this, part, theData->Table(i), dir, tnames[i],
1147 verboseLevel, ascii, spline)) {
1148 return false;
1149 }
1150 }
1151 return true;
1152}
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 * G4VEnergyLossProcess::SecondaryParticle ( ) const
inline

Definition at line 860 of file G4VEnergyLossProcess.hh.

861{
862 return secondaryParticle;
863}

◆ SelectModel()

void G4VEnergyLossProcess::SelectModel ( G4double kinEnergy)
inlineprotected

Definition at line 544 of file G4VEnergyLossProcess.hh.

545{
546 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
547 currentModel->SetCurrentCouple(currentCouple);
548}

Referenced by CrossSectionPerVolume(), GetDEDXDispersion(), PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ SelectModelForMaterial()

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

Definition at line 552 of file G4VEnergyLossProcess.hh.

554{
555 return modelManager->SelectModel(kinEnergy, idx);
556}

◆ SetBaseParticle()

void G4VEnergyLossProcess::SetBaseParticle ( const G4ParticleDefinition * p)
inline

Definition at line 838 of file G4VEnergyLossProcess.hh.

839{
840 baseParticle = p;
841}

Referenced by G4hIonisation::InitialiseEnergyLossProcess(), and G4ionIonisation::InitialiseEnergyLossProcess().

◆ SetCrossSectionBiasingFactor()

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

Definition at line 1350 of file G4VEnergyLossProcess.cc.

1352{
1353 if(f > 0.0) {
1354 biasFactor = f;
1355 weightFlag = flag;
1356 if(1 < verboseLevel) {
1357 G4cout << "### SetCrossSectionBiasingFactor: for "
1358 << " process " << GetProcessName()
1359 << " biasFactor= " << f << " weightFlag= " << flag
1360 << G4endl;
1361 }
1362 }
1363}

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ SetCrossSectionType()

void G4VEnergyLossProcess::SetCrossSectionType ( G4CrossSectionType val)
inline

Definition at line 874 of file G4VEnergyLossProcess.hh.

875{
876 fXSType = val;
877}

Referenced by G4EmTableUtil::BuildLocalElossProcess(), and G4eBremsstrahlung::G4eBremsstrahlung().

◆ SetCSDARangeTable()

void G4VEnergyLossProcess::SetCSDARangeTable ( G4PhysicsTable * pRange)

Definition at line 1276 of file G4VEnergyLossProcess.cc.

1277{
1278 theCSDARangeTable = p;
1279}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetDEDXBinning()

void G4VEnergyLossProcess::SetDEDXBinning ( G4int nbins)

Definition at line 1458 of file G4VEnergyLossProcess.cc.

1459{
1460 if(2 < n && n < 1000000000) {
1461 nBins = n;
1462 actBinning = true;
1463 } else {
1464 G4double e = (G4double)n;
1465 PrintWarning("SetDEDXBinning", e);
1466 }
1467}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), and G4mplIonisation::InitialiseEnergyLossProcess().

◆ SetDEDXTable()

void G4VEnergyLossProcess::SetDEDXTable ( G4PhysicsTable * p,
G4EmTableType tType )

Definition at line 1250 of file G4VEnergyLossProcess.cc.

1251{
1252 if(1 < verboseLevel) {
1253 G4cout << "### Set DEDX table " << p << " " << theDEDXTable
1254 << " " << theDEDXunRestrictedTable << " " << theIonisationTable
1255 << " for " << particle->GetParticleName()
1256 << " and process " << GetProcessName()
1257 << " type=" << tType << " isIonisation:" << isIonisation << G4endl;
1258 }
1259 if(fTotal == tType) {
1260 theDEDXunRestrictedTable = p;
1261 } else if(fRestricted == tType) {
1262 theDEDXTable = p;
1263 if(isMaster && nullptr == baseParticle) {
1264 theData->UpdateTable(theDEDXTable, 0);
1265 }
1266 } else if(fIsIonisation == tType) {
1267 theIonisationTable = p;
1268 if(isMaster && nullptr == baseParticle) {
1269 theData->UpdateTable(theIonisationTable, 1);
1270 }
1271 }
1272}
@ fIsIonisation

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetDynamicMassCharge()

void G4VEnergyLossProcess::SetDynamicMassCharge ( G4double massratio,
G4double charge2ratio )

Definition at line 190 of file G4VEnergyLossProcess.cc.

192{
193 massRatio = massratio;
194 logMassRatio = G4Log(massRatio);
195 fFactor = charge2ratio*biasFactor;
196 if(baseMat) { fFactor *= (*theDensityFactor)[currentCoupleIndex]; }
197 chargeSqRatio = charge2ratio;
198 reduceFactor = 1.0/(fFactor*massRatio);
199}

◆ SetEmModel()

◆ SetEnergyOfCrossSectionMax()

void G4VEnergyLossProcess::SetEnergyOfCrossSectionMax ( std::vector< G4double > * p)

Definition at line 1328 of file G4VEnergyLossProcess.cc.

1329{
1330 theEnergyOfCrossSectionMax = p;
1331}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetFluctModel()

◆ SetInverseRangeTable()

void G4VEnergyLossProcess::SetInverseRangeTable ( G4PhysicsTable * p)

Definition at line 1290 of file G4VEnergyLossProcess.cc.

1291{
1292 theInverseRangeTable = p;
1293}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetIonisation()

◆ SetLambdaTable()

void G4VEnergyLossProcess::SetLambdaTable ( G4PhysicsTable * p)

Definition at line 1297 of file G4VEnergyLossProcess.cc.

1298{
1299 if(1 < verboseLevel) {
1300 G4cout << "### Set Lambda table " << p << " " << theLambdaTable
1301 << " for " << particle->GetParticleName()
1302 << " and process " << GetProcessName() << G4endl;
1303 }
1304 theLambdaTable = p;
1305 tablesAreBuilt = true;
1306
1307 if(isMaster && nullptr != p) {
1308 delete theEnergyOfCrossSectionMax;
1309 theEnergyOfCrossSectionMax = nullptr;
1310 if(fEmTwoPeaks == fXSType) {
1311 if(nullptr != fXSpeaks) {
1312 for(auto & ptr : *fXSpeaks) { delete ptr; }
1313 delete fXSpeaks;
1314 }
1315 G4LossTableBuilder* bld = lManager->GetTableBuilder();
1316 fXSpeaks = G4EmUtility::FillPeaksStructure(p, bld);
1317 if(nullptr == fXSpeaks) { fXSType = fEmOnePeak; }
1318 }
1319 if(fXSType == fEmOnePeak) {
1320 theEnergyOfCrossSectionMax = G4EmUtility::FindCrossSectionMax(p);
1321 if(nullptr == theEnergyOfCrossSectionMax) { fXSType = fEmIncreasing; }
1322 }
1323 }
1324}
@ fEmOnePeak
@ fEmTwoPeaks
@ fEmIncreasing
static std::vector< G4TwoPeaksXS * > * FillPeaksStructure(G4PhysicsTable *, G4LossTableBuilder *)
static std::vector< G4double > * FindCrossSectionMax(G4PhysicsTable *)

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetLinearLossLimit()

void G4VEnergyLossProcess::SetLinearLossLimit ( G4double val)

Definition at line 1427 of file G4VEnergyLossProcess.cc.

1428{
1429 if(0.0 < val && val < 1.0) {
1430 linLossLimit = val;
1431 actLinLossLimit = true;
1432 } else { PrintWarning("SetLinearLossLimit", val); }
1433}

Referenced by G4ionIonisation::G4ionIonisation().

◆ SetLossFluctuations()

void G4VEnergyLossProcess::SetLossFluctuations ( G4bool )

Definition at line 1419 of file G4VEnergyLossProcess.cc.

1420{
1421 G4cout << "### G4VEnergyLossProcess::SetLossFluctuations has no effect and "
1422 << "will be removed for the next major release" << G4endl;
1423}

◆ SetLowestEnergyLimit()

void G4VEnergyLossProcess::SetLowestEnergyLimit ( G4double val)

Definition at line 1450 of file G4VEnergyLossProcess.cc.

1451{
1452 if(1.e-18 < val && val < 1.e+50) { lowestKinEnergy = val; }
1453 else { PrintWarning("SetLowestEnergyLimit", val); }
1454}

◆ SetMaxKinEnergy()

void G4VEnergyLossProcess::SetMaxKinEnergy ( G4double e)

Definition at line 1481 of file G4VEnergyLossProcess.cc.

1482{
1483 if(minKinEnergy < e && e < 1.e+50) {
1484 maxKinEnergy = e;
1485 actMaxKinEnergy = true;
1486 if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
1487 } else { PrintWarning("SetMaxKinEnergy", e); }
1488}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), and G4mplIonisation::InitialiseEnergyLossProcess().

◆ SetMinKinEnergy()

void G4VEnergyLossProcess::SetMinKinEnergy ( G4double e)

Definition at line 1471 of file G4VEnergyLossProcess.cc.

1472{
1473 if(1.e-18 < e && e < maxKinEnergy) {
1474 minKinEnergy = e;
1475 actMinKinEnergy = true;
1476 } else { PrintWarning("SetMinKinEnergy", e); }
1477}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), and G4mplIonisation::InitialiseEnergyLossProcess().

◆ SetParticle()

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

Definition at line 822 of file G4VEnergyLossProcess.hh.

823{
824 particle = p;
825}

◆ SetRangeTableForLoss()

void G4VEnergyLossProcess::SetRangeTableForLoss ( G4PhysicsTable * p)

Definition at line 1283 of file G4VEnergyLossProcess.cc.

1284{
1285 theRangeTableForLoss = p;
1286}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ SetSecondaryParticle()

◆ SetSpline()

void G4VEnergyLossProcess::SetSpline ( G4bool val)
inline

Definition at line 867 of file G4VEnergyLossProcess.hh.

868{
869 spline = val;
870}

Referenced by G4ePairProduction::G4ePairProduction().

◆ SetStepFunction()

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

Definition at line 1437 of file G4VEnergyLossProcess.cc.

1438{
1439 if(0.0 < v1 && 1.0 >= v1 && 0.0 < v2) {
1440 dRoverRange = std::min(1.0, v1);
1441 finalRange = std::min(v2, 1.e+50);
1442 } else {
1443 PrintWarning("SetStepFunctionV1", v1);
1444 PrintWarning("SetStepFunctionV2", v2);
1445 }
1446}

Referenced by LBE::ConstructEM(), and G4EmExtraParameters::FillStepFunction().

◆ SetTwoPeaksXS()

void G4VEnergyLossProcess::SetTwoPeaksXS ( std::vector< G4TwoPeaksXS * > * ptr)

Definition at line 1335 of file G4VEnergyLossProcess.cc.

1336{
1337 fXSpeaks = ptr;
1338}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ StartTracking()

void G4VEnergyLossProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 535 of file G4VEnergyLossProcess.cc.

536{
537 // reset parameters for the new track
540 preStepLambda = 0.0;
541 currentCouple = nullptr;
542
543 // reset ion
544 if(isIon) {
545 const G4double newmass = track->GetParticleDefinition()->GetPDGMass();
546 massRatio = (nullptr == baseParticle) ? CLHEP::proton_mass_c2/newmass
547 : baseParticle->GetPDGMass()/newmass;
548 logMassRatio = G4Log(massRatio);
549 }
550 // forced biasing only for primary particles
551 if(nullptr != biasManager) {
552 if(0 == track->GetParentID()) {
553 biasFlag = true;
554 biasManager->ResetForcedInteraction();
555 }
556 }
557}

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 1113 of file G4VEnergyLossProcess.cc.

1115{
1116 if (!isMaster || nullptr != baseParticle || part != particle ) return true;
1117 for(std::size_t i=0; i<7; ++i) {
1118 // ionisation table only for ionisation process
1119 if (nullptr == theData->Table(i) || (!isIonisation && 1 == i)) {
1120 continue;
1121 }
1122 if (-1 < verboseLevel) {
1123 G4cout << "G4VEnergyLossProcess::StorePhysicsTable i=" << i
1124 << " " << particle->GetParticleName()
1125 << " " << GetProcessName()
1126 << " " << tnames[i] << " " << theData->Table(i) << G4endl;
1127 }
1128 if (!G4EmTableUtil::StoreTable(this, part, theData->Table(i),
1129 dir, tnames[i], verboseLevel, ascii)) {
1130 return false;
1131 }
1132 }
1133 return true;
1134}
static G4bool StoreTable(G4VProcess *, const G4ParticleDefinition *, G4PhysicsTable *, const G4String &dir, const G4String &tname, G4int verb, G4bool ascii)

◆ StreamProcessInfo()

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

Reimplemented in G4eBremsstrahlung, G4ePairProduction, and G4MuPairProduction.

Definition at line 98 of file G4VEnergyLossProcess.hh.

98{};

◆ TablesAreBuilt()

G4bool G4VEnergyLossProcess::TablesAreBuilt ( ) const
inline

Definition at line 923 of file G4VEnergyLossProcess.hh.

924{
925 return tablesAreBuilt;
926}

◆ TwoPeaksXS()

std::vector< G4TwoPeaksXS * > * G4VEnergyLossProcess::TwoPeaksXS ( ) const
inline

Definition at line 994 of file G4VEnergyLossProcess.hh.

995{
996 return fXSpeaks;
997}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

◆ UseBaseMaterial()

G4bool G4VEnergyLossProcess::UseBaseMaterial ( ) const
inline

Definition at line 979 of file G4VEnergyLossProcess.hh.

980{
981 return baseMat;
982}

Referenced by G4EmTableUtil::BuildLocalElossProcess().

Member Data Documentation

◆ currentCouple

const G4MaterialCutsCouple* G4VEnergyLossProcess::currentCouple = nullptr
protected

◆ currentCoupleIndex

◆ currentMaterial

const G4Material* G4VEnergyLossProcess::currentMaterial = nullptr
protected

◆ fParticleChange

G4ParticleChangeForLoss G4VEnergyLossProcess::fParticleChange
protected

Definition at line 418 of file G4VEnergyLossProcess.hh.

Referenced by AlongStepDoIt(), G4VEnergyLossProcess(), and PostStepDoIt().

◆ mfpKinEnergy

G4double G4VEnergyLossProcess::mfpKinEnergy = 0.0
protected

◆ preStepKinEnergy

G4double G4VEnergyLossProcess::preStepKinEnergy = 0.0
protected

Definition at line 480 of file G4VEnergyLossProcess.hh.

Referenced by AlongStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ preStepLambda

G4double G4VEnergyLossProcess::preStepLambda = 0.0
protected

◆ preStepScaledEnergy

G4double G4VEnergyLossProcess::preStepScaledEnergy = 0.0
protected

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