Geant4  10.02.p03
G4hhIonisation Class Reference

#include <G4hhIonisation.hh>

Inheritance diagram for G4hhIonisation:
Collaboration diagram for G4hhIonisation:

Public Member Functions

 G4hhIonisation (const G4String &name="hhIoni")
 
virtual ~G4hhIonisation ()
 
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

G4hhIonisationoperator= (const G4hhIonisation &right)
 
 G4hhIonisation (const G4hhIonisation &)
 

Private Attributes

G4double mass
 
G4double ratio
 
const G4ParticleDefinitiontheParticle
 
const G4ParticleDefinitiontheBaseParticle
 
G4VEmFluctuationModelflucModel
 
G4bool isInitialised
 

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 61 of file G4hhIonisation.hh.

Constructor & Destructor Documentation

◆ G4hhIonisation() [1/2]

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

Definition at line 64 of file G4hhIonisation.cc.

65  : G4VEnergyLossProcess(name),
66  theParticle(0),
67  theBaseParticle(0),
68  isInitialised(false)
69 {
70  SetStepFunction(0.1, 0.1*mm);
71  SetVerboseLevel(1);
74  mass = 0.0;
75  ratio = 0.0;
76  flucModel = 0;
77 }
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4VEmFluctuationModel * flucModel
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
const G4ParticleDefinition * theBaseParticle
const G4ParticleDefinition * theParticle
static const double mm
Definition: G4SIunits.hh:114
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
Here is the call graph for this function:

◆ ~G4hhIonisation()

G4hhIonisation::~G4hhIonisation ( )
virtual

Definition at line 81 of file G4hhIonisation.cc.

82 {}

◆ G4hhIonisation() [2/2]

G4hhIonisation::G4hhIonisation ( const G4hhIonisation )
private

Member Function Documentation

◆ InitialiseEnergyLossProcess()

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

Implements G4VEnergyLossProcess.

Definition at line 106 of file G4hhIonisation.cc.

109 {
110  if(isInitialised) { return; }
111 
112  theParticle = part;
113  if(bpart) {
114  G4cout << "G4hhIonisation::InitialiseEnergyLossProcess WARNING: no "
115  << "base particle should be defined for the process "
116  << GetProcessName() << G4endl;
117  }
118  SetBaseParticle(0);
123 
125  G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
126  G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
127 
128  SetMinKinEnergy(emin);
129  SetMaxKinEnergy(emax);
130  G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
131  SetDEDXBinning(bin);
132 
133  G4VEmModel* em = 0;
134  if(part->GetPDGCharge() > 0.0) { em = new G4BraggNoDeltaModel(); }
135  else { em = new G4ICRU73NoDeltaModel(); }
136  em->SetLowEnergyLimit(emin);
137  em->SetHighEnergyLimit(eth);
138  AddEmModel(1, em, flucModel);
139 
140  em = new G4BetheBlochNoDeltaModel();
141  em->SetLowEnergyLimit(eth);
142  em->SetHighEnergyLimit(emax);
143  AddEmModel(1, em, flucModel);
144 
145  if(verboseLevel>1) {
146  G4cout << "G4hhIonisation is initialised" << G4endl;
147  }
148  isInitialised = true;
149 }
static const double MeV
Definition: G4SIunits.hh:211
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VEmFluctuationModel * flucModel
float bin[41]
Definition: plottest35.C:14
G4int NumberOfBinsPerDecade() const
int G4int
Definition: G4Types.hh:78
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4double MinKinEnergy() const
TString part[npart]
float proton_mass_c2
Definition: hepunit.py:275
void SetMaxKinEnergy(G4double e)
float electron_mass_c2
Definition: hepunit.py:274
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
static const G4double emax
G4double MaxKinEnergy() const
int G4lrint(double ad)
Definition: templates.hh:163
static G4EmParameters * Instance()
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * theParticle
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4double GetPDGCharge() const
void SetBaseParticle(const G4ParticleDefinition *p)
void SetDEDXBinning(G4int nbins)
void SetMinKinEnergy(G4double e)
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4hhIonisation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEnergyLossProcess.

Definition at line 86 of file G4hhIonisation.cc.

87 {
88  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 100.0*MeV &&
89  !p.IsShortLived());
90 }
static const double MeV
Definition: G4SIunits.hh:211
G4double GetPDGCharge() const
Here is the call graph for this function:

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 94 of file G4hhIonisation.cc.

97 {
98  G4double x = 0.5*cut/electron_mass_c2;
100  G4double gam = x*y + std::sqrt((1. + x)*(1. + x*y*y));
101  return mass*(gam - 1.0);
102 }
Double_t y
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76

◆ operator=()

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

◆ PrintInfo()

void G4hhIonisation::PrintInfo ( )
virtual

Implements G4VEnergyLossProcess.

Definition at line 153 of file G4hhIonisation.cc.

154 {
155  G4cout << " Delta-ray will not be produced; "
156  << G4endl;
157 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Data Documentation

◆ flucModel

G4VEmFluctuationModel* G4hhIonisation::flucModel
private

Definition at line 94 of file G4hhIonisation.hh.

◆ isInitialised

G4bool G4hhIonisation::isInitialised
private

Definition at line 96 of file G4hhIonisation.hh.

◆ mass

G4double G4hhIonisation::mass
private

Definition at line 89 of file G4hhIonisation.hh.

◆ ratio

G4double G4hhIonisation::ratio
private

Definition at line 90 of file G4hhIonisation.hh.

◆ theBaseParticle

const G4ParticleDefinition* G4hhIonisation::theBaseParticle
private

Definition at line 93 of file G4hhIonisation.hh.

◆ theParticle

const G4ParticleDefinition* G4hhIonisation::theParticle
private

Definition at line 92 of file G4hhIonisation.hh.


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