Geant4 11.1.1
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")
 
 ~G4hIonisation () override
 
G4bool IsApplicable (const G4ParticleDefinition &p) override
 
G4double MinPrimaryEnergy (const G4ParticleDefinition *p, const G4Material *, G4double cut) final
 
virtual void ProcessDescription (std::ostream &) const override
 
G4hIonisationoperator= (const G4hIonisation &right)=delete
 
 G4hIonisation (const G4hIonisation &)=delete
 
- Public Member Functions inherited from G4VEnergyLossProcess
 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 val)
 
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 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 const G4VProcessGetCreatorProcess () const
 
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

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
 
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
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 = nullptr
 
const G4MaterialCutsCouplecurrentCouple = nullptr
 
G4double preStepLambda = 0.0
 
G4double preStepKinEnergy = 0.0
 
G4double preStepLogKinEnergy = LOG_EKIN_MIN
 
G4double preStepScaledEnergy = 0.0
 
G4double preStepLogScaledEnergy = LOG_EKIN_MIN
 
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
 

Detailed Description

Definition at line 84 of file G4hIonisation.hh.

Constructor & Destructor Documentation

◆ G4hIonisation() [1/2]

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

Definition at line 62 of file G4hIonisation.cc.

64 isInitialised(false)
65{
68 mass = 0.0;
69 ratio = 0.0;
70 eth = 2*CLHEP::MeV;
71}
@ fIonisation
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4hIonisation()

G4hIonisation::~G4hIonisation ( )
overridedefault

◆ G4hIonisation() [2/2]

G4hIonisation::G4hIonisation ( const G4hIonisation )
delete

Member Function Documentation

◆ InitialiseEnergyLossProcess()

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

Implements G4VEnergyLossProcess.

Definition at line 97 of file G4hIonisation.cc.

100{
101 if(!isInitialised) {
102
103 const G4ParticleDefinition* theBaseParticle = nullptr;
104 G4String pname = part->GetParticleName();
105 G4double q = part->GetPDGCharge();
106
107 //G4cout << " G4hIonisation::InitialiseEnergyLossProcess " << pname
108 // << " " << bpart << G4endl;
109
110 // define base particle
111 if(part == bpart) {
112 theBaseParticle = nullptr;
113 } else if(nullptr != bpart) {
114 theBaseParticle = bpart;
115
116 } else if(pname == "proton" || pname == "anti_proton" ||
117 pname == "pi+" || pname == "pi-" ||
118 pname == "kaon+" || pname == "kaon-" ||
119 pname == "GenericIon" || pname == "alpha") {
120 // no base particles
121 theBaseParticle = nullptr;
122
123 } else {
124 // select base particle
125 if(part->GetPDGSpin() == 0.0) {
126 if(q > 0.0) { theBaseParticle = G4KaonPlus::KaonPlus(); }
127 else { theBaseParticle = G4KaonMinus::KaonMinus(); }
128 } else {
129 if(q > 0.0) { theBaseParticle = G4Proton::Proton(); }
130 else { theBaseParticle = G4AntiProton::AntiProton(); }
131 }
132 }
133 SetBaseParticle(theBaseParticle);
134
135 // model limit defined for protons
136 mass = part->GetPDGMass();
137 ratio = electron_mass_c2/mass;
138 eth = 2.0*MeV*mass/proton_mass_c2;
139
141 G4double emin = param->MinKinEnergy();
142 G4double emax = param->MaxKinEnergy();
143
144 // define model of energy loss fluctuations
145 if (nullptr == FluctModel()) {
146 G4bool ion = (pname == "GenericIon" || pname == "alpha");
148 }
149
150 if (nullptr == EmModel(0)) {
151 if(q > 0.0) { SetEmModel(new G4BraggModel()); }
152 else { SetEmModel(new G4ICRU73QOModel()); }
153 }
154 // to compute ranges correctly we have to use low-energy
155 // model even if activation limit is high
156 EmModel(0)->SetLowEnergyLimit(emin);
157
158 // high energy limit may be eth or DBL_MAX
159 G4double emax1 = (EmModel(0)->HighEnergyLimit() < emax) ? eth : emax;
160 EmModel(0)->SetHighEnergyLimit(emax1);
161 AddEmModel(1, EmModel(0), FluctModel());
162
163 // second model is used if the first does not cover energy range
164 if(emax1 < emax) {
165 if (nullptr == EmModel(1)) { SetEmModel(new G4BetheBlochModel()); }
166 EmModel(1)->SetLowEnergyLimit(emax1);
167
168 // for extremely heavy particles upper limit of the model
169 // should be increased
170 emax = std::max(emax, eth*10);
171 EmModel(1)->SetHighEnergyLimit(emax);
172 AddEmModel(2, EmModel(1), FluctModel());
173 }
174 isInitialised = true;
175 }
176}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
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:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
G4VEmModel * EmModel(std::size_t index=0) const
void SetEmModel(G4VEmModel *, G4int index=0)
void SetBaseParticle(const G4ParticleDefinition *p)
G4VEmFluctuationModel * FluctModel() const

◆ IsApplicable()

G4bool G4hIonisation::IsApplicable ( const G4ParticleDefinition p)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 79 of file G4hIonisation.cc.

80{
81 return true;
82}

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 86 of file G4hIonisation.cc.

89{
90 G4double x = 0.5*cut/electron_mass_c2;
91 G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
92 return mass*(gam - 1.0);
93}

◆ operator=()

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

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 180 of file G4hIonisation.cc.

181{
182 out << " Ionisation";
184}
void ProcessDescription(std::ostream &outFile) const override

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