Geant4  10.02.p03
G4hIonisation Class Reference

#include <G4hIonisation.hh>

Inheritance diagram for G4hIonisation:
Collaboration 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 void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
void PrintInfoDefinition (const G4ParticleDefinition &part)
 
virtual void StartTracking (G4Track *)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
G4double SampleSubCutSecondaries (std::vector< G4Track *> &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual 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) 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=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 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 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
 
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 ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (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

virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
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
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

G4hIonisationoperator= (const G4hIonisation &right)
 
 G4hIonisation (const G4hIonisation &)
 

Private Attributes

G4bool isInitialised
 
G4double mass
 
G4double ratio
 
G4double eth
 

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
 
G4VParticleChange * pParticleChange
 
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() [1/2]

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

Definition at line 110 of file G4hIonisation.cc.

111  : G4VEnergyLossProcess(name),
112  isInitialised(false)
113 {
114  SetStepFunction(0.2, 0.1*mm);
117  mass = 0.0;
118  ratio = 0.0;
119  eth = 2*MeV;
120 }
G4bool isInitialised
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
static const double MeV
Definition: G4SIunits.hh:211
void SetStepFunction(G4double v1, G4double v2)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetSecondaryParticle(const G4ParticleDefinition *p)
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ ~G4hIonisation()

G4hIonisation::~G4hIonisation ( )
virtual

Definition at line 124 of file G4hIonisation.cc.

125 {}

◆ G4hIonisation() [2/2]

G4hIonisation::G4hIonisation ( const G4hIonisation )
private

Member Function Documentation

◆ InitialiseEnergyLossProcess()

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

Implements G4VEnergyLossProcess.

Definition at line 148 of file G4hIonisation.cc.

151 {
152  if(!isInitialised) {
153 
154  const G4ParticleDefinition* theBaseParticle = 0;
155  G4String pname = part->GetParticleName();
156  G4double q = part->GetPDGCharge();
157 
158  //G4cout << " G4hIonisation::InitialiseEnergyLossProcess " << pname
159  // << " " << bpart << G4endl;
160 
161  // standard base particles
162  if(part == bpart || pname == "proton" ||
163  pname == "anti_proton" ||
164  pname == "pi+" || pname == "pi-" ||
165  pname == "kaon+" || pname == "kaon-" || pname == "GenericIon"
166  || pname == "He3" || pname == "alpha")
167  {
168  theBaseParticle = 0;
169  }
170  // select base particle
171  else if(bpart == 0) {
172 
173  if(part->GetPDGSpin() == 0.0) {
174  if(q > 0.0) { theBaseParticle = G4KaonPlus::KaonPlus(); }
175  else { theBaseParticle = G4KaonMinus::KaonMinus(); }
176  } else {
177  if(q > 0.0) { theBaseParticle = G4Proton::Proton(); }
178  else { theBaseParticle = G4AntiProton::AntiProton(); }
179  }
180 
181  // base particle defined by interface
182  } else {
183  theBaseParticle = bpart;
184  }
185  SetBaseParticle(theBaseParticle);
186 
187  mass = part->GetPDGMass();
189  eth = 2.0*MeV*mass/proton_mass_c2;
190 
191  if (!EmModel(1)) {
192  if(q > 0.0) { SetEmModel(new G4BraggModel(),1); }
193  else { SetEmModel(new G4ICRU73QOModel(),1); }
194  }
196  EmModel(1)->SetLowEnergyLimit(param->MinKinEnergy());
198  AddEmModel(1, EmModel(1), new G4IonFluctuations());
199 
201 
202  if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
205  AddEmModel(2, EmModel(2), FluctModel());
206 
207  isInitialised = true;
208  }
209 }
G4bool isInitialised
static const double MeV
Definition: G4SIunits.hh:211
void SetFluctModel(G4VEmFluctuationModel *)
G4VEmModel * EmModel(G4int index=1) const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4VEmFluctuationModel * FluctModel()
const G4String & GetParticleName() const
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
G4double MinKinEnergy() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
string pname
Definition: eplot.py:33
G4double MaxKinEnergy() const
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double GetPDGCharge() const
void SetBaseParticle(const G4ParticleDefinition *p)
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4hIonisation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEnergyLossProcess.

Definition at line 129 of file G4hIonisation.cc.

130 {
131  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV &&
132  !p.IsShortLived());
133 }
static const double MeV
Definition: G4SIunits.hh:211
G4double GetPDGCharge() const
Here is the call graph for this function:

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 137 of file G4hIonisation.cc.

140 {
141  G4double x = 0.5*cut/electron_mass_c2;
142  G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
143  return mass*(gam - 1.0);
144 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76

◆ operator=()

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

◆ PrintInfo()

void G4hIonisation::PrintInfo ( )
virtual

Implements G4VEnergyLossProcess.

Definition at line 213 of file G4hIonisation.cc.

214 {}

Member Data Documentation

◆ eth

G4double G4hIonisation::eth
private

Definition at line 117 of file G4hIonisation.hh.

◆ isInitialised

G4bool G4hIonisation::isInitialised
private

Definition at line 113 of file G4hIonisation.hh.

◆ mass

G4double G4hIonisation::mass
private

Definition at line 115 of file G4hIonisation.hh.

◆ ratio

G4double G4hIonisation::ratio
private

Definition at line 116 of file G4hIonisation.hh.


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