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

#include <G4NuclearStopping.hh>

+ Inheritance diagram for G4NuclearStopping:

Public Member Functions

 G4NuclearStopping (const G4String &processName="nuclearStopping")
 
virtual ~G4NuclearStopping ()
 
G4bool IsApplicable (const G4ParticleDefinition &p)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step)
 
void PrintInfo ()
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double &kinEnergy, const G4MaterialCutsCouple *couple)
 
void SetLambdaBinning (G4int nbins)
 
G4int LambdaBinning () const
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMinKinEnergyPrim (G4double e)
 
const G4PhysicsTableLambdaTable () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
void SetModel (G4VEmModel *, G4int index=1)
 
G4VEmModelModel (G4int index=1)
 
G4VEmModelEmModel (G4int index=1)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
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 SetPolarAngleLimit (G4double a)
 
G4double PolarAngleLimit () const
 
void SetLambdaFactor (G4double val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetApplyCuts (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
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 void StartTracking (G4Track *)
 
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
 

Protected Member Functions

void InitialiseProcess (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEmProcess
G4ParticleChangeForGamma fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 65 of file G4NuclearStopping.hh.

Constructor & Destructor Documentation

◆ G4NuclearStopping()

G4NuclearStopping::G4NuclearStopping ( const G4String processName = "nuclearStopping")

Definition at line 54 of file G4NuclearStopping.cc.

55 : G4VEmProcess(processName)
56{
57 isInitialized = false;
59 SetBuildTableFlag(false);
60 enableAlongStepDoIt = true;
61 enablePostStepDoIt = false;
62 modelICRU49 = 0;
63}
@ fNuclearStopping
void SetBuildTableFlag(G4bool val)
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352

◆ ~G4NuclearStopping()

G4NuclearStopping::~G4NuclearStopping ( )
virtual

Definition at line 67 of file G4NuclearStopping.cc.

68{}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4NuclearStopping::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 105 of file G4NuclearStopping.cc.

107{
108 nParticleChange.InitializeForAlongStep(track);
110
111 const G4ParticleDefinition* part = track.GetParticleDefinition();
112 G4double Z = std::fabs(part->GetPDGCharge()/eplus);
113
114 if(T2 > 0.0 && T2*proton_mass_c2 < Z*Z*MeV*part->GetPDGMass()) {
115
116 G4double length = step.GetStepLength();
117 if(length > 0.0) {
118
119 // primary
121 G4double T = 0.5*(T1 + T2);
122 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
123 G4VEmModel* mod = SelectModel(T, couple->GetIndex());
124
125 // sample stopping
126 if(modelICRU49) { modelICRU49->SetFluctuationFlag(true); }
127 G4double nloss =
128 length*mod->ComputeDEDXPerVolume(couple->GetMaterial(), part, T);
129 if(nloss > T1) { nloss = T1; }
130 nParticleChange.SetProposedKineticEnergy(T1 - nloss);
131 nParticleChange.ProposeLocalEnergyDeposit(nloss);
132 nParticleChange.ProposeNonIonizingEnergyDeposit(nloss);
133 if(modelICRU49) { modelICRU49->SetFluctuationFlag(false); }
134 }
135 }
136 return &nParticleChange;
137}
double G4double
Definition: G4Types.hh:64
const G4Material * GetMaterial() const
void InitializeForAlongStep(const G4Track &)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double GetPDGCharge() const
G4double GetKineticEnergy() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:177
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4NuclearStopping::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 96 of file G4NuclearStopping.cc.

98{
99 *selection = NotCandidateForSelection;
100 return DBL_MAX;
101}
@ NotCandidateForSelection
#define DBL_MAX
Definition: templates.hh:83

◆ InitialiseProcess()

void G4NuclearStopping::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 79 of file G4NuclearStopping.cc.

80{
81 if(!isInitialized) {
82 isInitialized = true;
83
84 if(!EmModel(1)) {
85 modelICRU49 = new G4ICRU49NuclearStoppingModel();
86 SetEmModel(modelICRU49);
87 }
88 AddEmModel(1, EmModel());
89 EmModel()->SetParticleChange(&nParticleChange);
90 }
91}
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:318
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)

◆ IsApplicable()

G4bool G4NuclearStopping::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 72 of file G4NuclearStopping.cc.

73{
74 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
75}

◆ PrintInfo()

void G4NuclearStopping::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 141 of file G4NuclearStopping.cc.

142{}

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