Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NuclearStopping Class Reference

#include <G4NuclearStopping.hh>

Inheritance diagram for G4NuclearStopping:
Collaboration diagram for G4NuclearStopping:

Public Member Functions

 G4NuclearStopping (const G4String &processName="nuclearStopping")
 
virtual ~G4NuclearStopping ()
 
G4bool IsApplicable (const G4ParticleDefinition &p) final
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) final
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step) final
 
virtual void PrintInfo () final
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
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)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4int GetNumberOfModels () const
 
G4int GetNumberOfRegionModels (size_t couple_index) const
 
G4VEmModelGetRegionModel (G4int idx, size_t couple_index) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4VEmModelGetCurrentModel () 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 SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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)
 
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 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

void InitialiseProcess (const G4ParticleDefinition *) final
 
- Protected Member Functions inherited from G4VEmProcess
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
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)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
- 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 64 of file G4NuclearStopping.hh.

Constructor & Destructor Documentation

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

Definition at line 56 of file G4NuclearStopping.cc.

57  : G4VEmProcess(processName)
58 {
59  isInitialized = false;
61  SetBuildTableFlag(false);
62  enableAlongStepDoIt = true;
63  enablePostStepDoIt = false;
64 }
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:91
void SetBuildTableFlag(G4bool val)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351

Here is the call graph for this function:

G4NuclearStopping::~G4NuclearStopping ( )
virtual

Definition at line 68 of file G4NuclearStopping.cc.

69 {}

Member Function Documentation

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

Reimplemented from G4VDiscreteProcess.

Definition at line 105 of file G4NuclearStopping.cc.

107 {
108  nParticleChange.InitializeForAlongStep(track);
109 
110  // this line only valid if nuclear stopping
111  // is computed after G4ionIonisation process
112  nParticleChange.SetProposedCharge(step.GetPostStepPoint()->GetCharge());
113 
115 
116  const G4ParticleDefinition* part = track.GetParticleDefinition();
117  G4double Z = std::abs(part->GetPDGCharge()/eplus);
118 
119  if(T2 > 0.0 && T2*proton_mass_c2 < Z*Z*MeV*part->GetPDGMass()) {
120 
121  G4double length = step.GetStepLength();
122  if(length > 0.0) {
123 
124  // primary
126  G4double T = 0.5*(T1 + T2);
127  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
128  G4VEmModel* mod = SelectModel(T, couple->GetIndex());
129 
130  // sample stopping
131  G4double nloss =
132  length*mod->ComputeDEDXPerVolume(couple->GetMaterial(), part, T);
133  if(nloss > T1) { nloss = T1; }
134  nParticleChange.SetProposedKineticEnergy(T1 - nloss);
135  nParticleChange.ProposeLocalEnergyDeposit(nloss);
136  nParticleChange.ProposeNonIonizingEnergyDeposit(nloss);
137  }
138  }
139  return &nParticleChange;
140 }
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4double GetStepLength() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void InitializeForAlongStep(const G4Track &)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4StepPoint * GetPreStepPoint() const
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
static constexpr double eplus
Definition: G4SIunits.hh:199
G4double GetCharge() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4ParticleDefinition * GetParticleDefinition() const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:248
G4StepPoint * GetPostStepPoint() const
void SetProposedCharge(G4double theCharge)
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 96 of file G4NuclearStopping.cc.

98 {
99  *selection = NotCandidateForSelection;
100  return DBL_MAX;
101 }
#define DBL_MAX
Definition: templates.hh:83
void G4NuclearStopping::InitialiseProcess ( const G4ParticleDefinition )
finalprotectedvirtual

Implements G4VEmProcess.

Definition at line 80 of file G4NuclearStopping.cc.

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

Here is the call graph for this function:

G4bool G4NuclearStopping::IsApplicable ( const G4ParticleDefinition p)
finalvirtual

Implements G4VEmProcess.

Definition at line 73 of file G4NuclearStopping.cc.

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

Here is the call graph for this function:

void G4NuclearStopping::PrintInfo ( )
finalvirtual

Implements G4VEmProcess.

Definition at line 144 of file G4NuclearStopping.cc.

145 {}

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