Geant4 9.6.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)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *p, const G4Material *, G4double cut)
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEnergyLossProcess
 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEnergyLossProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4double SampleSubCutSecondaries (std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
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 GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
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=0)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
G4VEmModelEmModel (G4int index=1)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
G4int NumberOfModels ()
 
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=0)
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &region="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void AddCollaborativeProcess (G4VEnergyLossProcess *)
 
void SetLossFluctuations (G4bool val)
 
void SetRandomStep (G4bool val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetIonisation (G4bool val)
 
G4bool IsIonisationProcess () const
 
void SetLinearLossLimit (G4double val)
 
void SetMinSubRange (G4double val)
 
void SetLambdaFactor (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 SetLambdaBinning (G4int nbins)
 
void SetDEDXBinningForCSDARange (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMaxKinEnergyForCSDARange (G4double e)
 
G4double CrossSectionBiasingFactor () const
 
G4double GetDEDX (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetDEDXForSubsec (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetCSDARange (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetKineticEnergy (G4double &range, const G4MaterialCutsCouple *)
 
G4double GetLambda (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4bool TablesAreBuilt () const
 
G4PhysicsTableDEDXTable () const
 
G4PhysicsTableDEDXTableForSubsec () const
 
G4PhysicsTableDEDXunRestrictedTable () const
 
G4PhysicsTableIonisationTable () const
 
G4PhysicsTableIonisationTableForSubsec () const
 
G4PhysicsTableCSDARangeTable () const
 
G4PhysicsTableRangeTableForLoss () const
 
G4PhysicsTableInverseRangeTable () const
 
G4PhysicsTableLambdaTable ()
 
G4PhysicsTableSubLambdaTable ()
 
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 ()
 
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 ()
 
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

virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetCurrentRange () 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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEnergyLossProcess
G4ParticleChangeForLoss 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 85 of file G4hIonisation.hh.

Constructor & Destructor Documentation

◆ G4hIonisation()

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

Definition at line 109 of file G4hIonisation.cc.

110 : G4VEnergyLossProcess(name),
111 isInitialised(false)
112{
113 SetStepFunction(0.2, 0.1*mm);
116 mass = 0.0;
117 ratio = 0.0;
118 eth = 2*MeV;
119}
@ fIonisation
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetStepFunction(G4double v1, G4double v2)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4hIonisation()

G4hIonisation::~G4hIonisation ( )
virtual

Definition at line 123 of file G4hIonisation.cc.

124{}

Member Function Documentation

◆ InitialiseEnergyLossProcess()

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

Implements G4VEnergyLossProcess.

Definition at line 147 of file G4hIonisation.cc.

150{
151 if(!isInitialised) {
152
153 const G4ParticleDefinition* theBaseParticle = 0;
154 G4String pname = part->GetParticleName();
155 G4double q = part->GetPDGCharge();
156
157 // standard base particles
158 if(part == bpart || pname == "proton" ||
159 pname == "anti_proton" ||
160 pname == "pi+" || pname == "pi-" ||
161 pname == "kaon+" || pname == "kaon-")
162 {
163 theBaseParticle = 0;
164 }
165 // select base particle
166 else if(bpart == 0) {
167
168 if(part->GetPDGSpin() == 0.0) {
169 if(q > 0.0) { theBaseParticle = G4KaonPlus::KaonPlus(); }
170 else { theBaseParticle = G4KaonMinus::KaonMinus(); }
171 } else {
172 if(q > 0.0) { theBaseParticle = G4Proton::Proton(); }
173 else { theBaseParticle = G4AntiProton::AntiProton(); }
174 }
175
176 // base particle defined by interface
177 } else {
178 theBaseParticle = bpart;
179 }
180 SetBaseParticle(theBaseParticle);
181
182 mass = part->GetPDGMass();
183 ratio = electron_mass_c2/mass;
184 eth = 2.0*MeV*mass/proton_mass_c2;
185
186 if (!EmModel(1)) {
187 if(q > 0.0) { SetEmModel(new G4BraggModel(),1); }
188 else { SetEmModel(new G4ICRU73QOModel(),1); }
189 }
191
192 // model limit defined for protons
193 //eth = (EmModel(1)->HighEnergyLimit())*mass/proton_mass_c2;
196
198
199 if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
200 EmModel(2)->SetLowEnergyLimit(eth);
202 AddEmModel(2, EmModel(2), FluctModel());
203
204 isInitialised = true;
205 }
206}
double G4double
Definition: G4Types.hh:64
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetBaseParticle(const G4ParticleDefinition *p)
G4VEmFluctuationModel * FluctModel()

◆ IsApplicable()

G4bool G4hIonisation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEnergyLossProcess.

Definition at line 128 of file G4hIonisation.cc.

129{
130 return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV &&
131 !p.IsShortLived());
132}

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 136 of file G4hIonisation.cc.

139{
140 G4double x = 0.5*cut/electron_mass_c2;
141 G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
142 return mass*(gam - 1.0);
143}

◆ PrintInfo()

void G4hIonisation::PrintInfo ( )
virtual

Implements G4VEnergyLossProcess.

Definition at line 210 of file G4hIonisation.cc.

211{}

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