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

#include <G4hhIonisation.hh>

+ Inheritance diagram for G4hhIonisation:

Public Member Functions

 G4hhIonisation (const G4String &name="hhIoni")
 
 ~G4hhIonisation () override
 
G4bool IsApplicable (const G4ParticleDefinition &p) override
 
G4double MinPrimaryEnergy (const G4ParticleDefinition *p, const G4Material *, G4double cut) override
 
void ProcessDescription (std::ostream &) const override
 
G4hhIonisationoperator= (const G4hhIonisation &right)=delete
 
 G4hhIonisation (const G4hhIonisation &)=delete
 
- Public Member Functions inherited from G4VEnergyLossProcess
 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 
 ~G4VEnergyLossProcess () override
 
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 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
 
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 const G4VProcessGetCreatorProcess () const
 
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 InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *) override
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual void StreamProcessInfo (std::ostream &) const
 
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
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 preStepScaledEnergy = 0.0
 
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 60 of file G4hhIonisation.hh.

Constructor & Destructor Documentation

◆ G4hhIonisation() [1/2]

G4hhIonisation::G4hhIonisation ( const G4String & name = "hhIoni")
explicit

Definition at line 62 of file G4hhIonisation.cc.

64{
68}
@ fIonisation
static G4Electron * Electron()
Definition G4Electron.cc:91
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)

◆ ~G4hhIonisation()

G4hhIonisation::~G4hhIonisation ( )
override

Definition at line 72 of file G4hhIonisation.cc.

73{}

◆ G4hhIonisation() [2/2]

G4hhIonisation::G4hhIonisation ( const G4hhIonisation & )
delete

Member Function Documentation

◆ InitialiseEnergyLossProcess()

void G4hhIonisation::InitialiseEnergyLossProcess ( const G4ParticleDefinition * part,
const G4ParticleDefinition * bpart )
overrideprotectedvirtual

Implements G4VEnergyLossProcess.

Definition at line 96 of file G4hhIonisation.cc.

99{
100 if(isInitialised) { return; }
101
102 theParticle = part;
103 if(nullptr != bpart) {
104 G4cout << "G4hhIonisation::InitialiseEnergyLossProcess WARNING: no "
105 << "base particle should be defined for the process "
106 << GetProcessName() << G4endl;
107 }
108 mass = theParticle->GetPDGMass();
109 ratio = electron_mass_c2/mass;
110 G4double eth = 2*CLHEP::MeV*mass/proton_mass_c2;
111 flucModel = new G4IonFluctuations();
112
114 G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
115 G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
116
117 SetMinKinEnergy(emin);
118 SetMaxKinEnergy(emax);
119 G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
120 SetDEDXBinning(bin);
121
122 G4VEmModel* em = EmModel(0);
123 if (nullptr == em) {
124 if(part->GetPDGCharge() > 0.0) { em = new G4BraggNoDeltaModel(); }
125 else { em = new G4ICRU73NoDeltaModel(); }
126 }
127 em->SetLowEnergyLimit(emin);
128 em->SetHighEnergyLimit(eth);
129 AddEmModel(1, em, flucModel);
130
131 em = EmModel(1);
132 if(nullptr == em) { em = new G4BetheBlochNoDeltaModel(); }
133 em->SetLowEnergyLimit(eth);
134 em->SetHighEnergyLimit(emax);
135 AddEmModel(1, em, flucModel);
136
137 if(verboseLevel>1) {
138 G4cout << "G4hhIonisation is initialised" << G4endl;
139 }
140 isInitialised = true;
141}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
G4VEmModel * EmModel(std::size_t index=0) const
void SetDEDXBinning(G4int nbins)
G4int verboseLevel
const G4String & GetProcessName() const
int G4lrint(double ad)
Definition templates.hh:134

◆ IsApplicable()

G4bool G4hhIonisation::IsApplicable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 77 of file G4hhIonisation.cc.

78{
79 return (p.GetPDGCharge() != 0.0);
80}

◆ MinPrimaryEnergy()

G4double G4hhIonisation::MinPrimaryEnergy ( const G4ParticleDefinition * p,
const G4Material * ,
G4double cut )
overridevirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 84 of file G4hhIonisation.cc.

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

◆ operator=()

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

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 145 of file G4hhIonisation.cc.

146{
147 out << "G4hhIonisation: no delta rays" << G4endl;
149}
void ProcessDescription(std::ostream &outFile) const override

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