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

#include <G4hIonisation.hh>

+ Inheritance diagram for G4hIonisation:

Public Member Functions

 G4hIonisation (const G4String &name="hIoni")
 
virtual ~G4hIonisation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *p, const G4Material *, G4double cut) final
 
virtual void PrintInfo () final
 
virtual void ProcessDescription (std::ostream &) const override
 
- Public Member Functions inherited from G4VEnergyLossProcess
 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEnergyLossProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
virtual void StartTracking (G4Track *) override
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double SampleSubCutSecondaries (std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
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 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, size_t &idx) const
 
void AddEmModel (G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4VEmModelEmModel (size_t index=0) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
G4int NumberOfModels () const
 
void SetFluctModel (G4VEmFluctuationModel *)
 
G4VEmFluctuationModelFluctModel ()
 
void SetBaseParticle (const G4ParticleDefinition *p)
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionBaseParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
void ActivateSubCutoff (G4bool val, const G4Region *region=nullptr)
 
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 AddCollaborativeProcess (G4VEnergyLossProcess *)
 
void SetLossFluctuations (G4bool val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () 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 SetSecondaryRangeTable (G4PhysicsTable *p)
 
void SetInverseRangeTable (G4PhysicsTable *p)
 
void SetLambdaTable (G4PhysicsTable *p)
 
void SetSubLambdaTable (G4PhysicsTable *p)
 
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 GetDEDX (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4double GetDEDXForSubsec (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4double GetCSDARange (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
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
 
G4PhysicsTableDEDXTableForSubsec () const
 
G4PhysicsTableDEDXunRestrictedTable () const
 
G4PhysicsTableIonisationTable () const
 
G4PhysicsTableIonisationTableForSubsec () const
 
G4PhysicsTableCSDARangeTable () const
 
G4PhysicsTableSecondaryRangeTable () const
 
G4PhysicsTableRangeTableForLoss () const
 
G4PhysicsTableInverseRangeTable () const
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableSubLambdaTable () 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 PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
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
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *) override
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
void SelectModel (G4double kinEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEnergyLossProcess
G4ParticleChangeForLoss fParticleChange
 
const G4MaterialcurrentMaterial
 
const G4MaterialCutsCouplecurrentCouple
 
size_t currentCoupleIndex
 
G4double preStepLambda
 
G4double fRange
 
G4double computedRange
 
G4double preStepKinEnergy
 
G4double preStepLogKinEnergy
 
G4double preStepScaledEnergy
 
G4double preStepLogScaledEnergy
 
G4double preStepRangeEnergy
 
G4double mfpKinEnergy
 
- 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
 

Detailed Description

Definition at line 84 of file G4hIonisation.hh.

Constructor & Destructor Documentation

◆ G4hIonisation()

G4hIonisation::G4hIonisation ( const G4String name = "hIoni")
explicit

Definition at line 69 of file G4hIonisation.cc.

71 isInitialised(false)
72{
75 mass = 0.0;
76 ratio = 0.0;
77 eth = 2*MeV;
78}
@ fIonisation
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4hIonisation()

G4hIonisation::~G4hIonisation ( )
virtual

Definition at line 82 of file G4hIonisation.cc.

83{}

Member Function Documentation

◆ InitialiseEnergyLossProcess()

void G4hIonisation::InitialiseEnergyLossProcess ( const G4ParticleDefinition part,
const G4ParticleDefinition bpart 
)
overrideprotectedvirtual

Implements G4VEnergyLossProcess.

Definition at line 106 of file G4hIonisation.cc.

109{
110 if(!isInitialised) {
111
112 const G4ParticleDefinition* theBaseParticle = nullptr;
113 G4String pname = part->GetParticleName();
114 G4double q = part->GetPDGCharge();
115
116 //G4cout << " G4hIonisation::InitialiseEnergyLossProcess " << pname
117 // << " " << bpart << G4endl;
118
119 // standard base particles
120 if(part == bpart || pname == "proton" ||
121 pname == "anti_proton" ||
122 pname == "pi+" || pname == "pi-" ||
123 pname == "kaon+" || pname == "kaon-" || pname == "GenericIon"
124 || pname == "He3" || pname == "alpha")
125 {
126 theBaseParticle = nullptr;
127 }
128 // select base particle
129 else if(bpart == nullptr) {
130
131 if(part->GetPDGSpin() == 0.0) {
132 if(q > 0.0) { theBaseParticle = G4KaonPlus::KaonPlus(); }
133 else { theBaseParticle = G4KaonMinus::KaonMinus(); }
134 } else {
135 if(q > 0.0) { theBaseParticle = G4Proton::Proton(); }
136 else { theBaseParticle = G4AntiProton::AntiProton(); }
137 }
138
139 // base particle defined by interface
140 } else {
141 theBaseParticle = bpart;
142 }
143 SetBaseParticle(theBaseParticle);
144
145 mass = part->GetPDGMass();
146 ratio = electron_mass_c2/mass;
147 eth = 2.0*MeV*mass/proton_mass_c2;
148
150 G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
151 G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
152
153 if(emin != param->MinKinEnergy() || emax != param->MaxKinEnergy()) {
154 SetMinKinEnergy(emin);
155 SetMaxKinEnergy(emax);
156 G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
157 SetDEDXBinning(bin);
158 }
159
160 if (!EmModel(0)) {
161 if(q > 0.0) { SetEmModel(new G4BraggModel()); }
162 else { SetEmModel(new G4ICRU73QOModel()); }
163 }
164 EmModel(0)->SetLowEnergyLimit(emin);
167
169
170 if (!EmModel(1)) { SetEmModel(new G4BetheBlochModel()); }
171 EmModel(1)->SetLowEnergyLimit(eth);
172 EmModel(1)->SetHighEnergyLimit(emax);
173 AddEmModel(1, EmModel(1), FluctModel());
174
175 isInitialised = true;
176 }
177}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
void SetDEDXBinning(G4int nbins)
void SetEmModel(G4VEmModel *, G4int index=0)
G4VEmModel * EmModel(size_t index=0) const
void SetBaseParticle(const G4ParticleDefinition *p)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4VEmFluctuationModel * FluctModel()
void SetMinKinEnergy(G4double e)
int G4lrint(double ad)
Definition: templates.hh:134

◆ IsApplicable()

G4bool G4hIonisation::IsApplicable ( const G4ParticleDefinition p)
overridevirtual

Implements G4VEnergyLossProcess.

Definition at line 87 of file G4hIonisation.cc.

88{
89 return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV &&
90 !p.IsShortLived());
91}

◆ MinPrimaryEnergy()

G4double G4hIonisation::MinPrimaryEnergy ( const G4ParticleDefinition p,
const G4Material ,
G4double  cut 
)
finalvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 95 of file G4hIonisation.cc.

98{
99 G4double x = 0.5*cut/electron_mass_c2;
100 G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
101 return mass*(gam - 1.0);
102}

◆ PrintInfo()

void G4hIonisation::PrintInfo ( )
finalvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 181 of file G4hIonisation.cc.

182{}

◆ ProcessDescription()

void G4hIonisation::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 186 of file G4hIonisation.cc.

187{
188 out << " Ionisation";
190}
virtual void ProcessDescription(std::ostream &outFile) const override

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