Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MuIonisation Class Reference

#include <G4MuIonisation.hh>

Inheritance diagram for G4MuIonisation:
Collaboration diagram for G4MuIonisation:

Public Member Functions

 G4MuIonisation (const G4String &name="muIoni")
 
virtual ~G4MuIonisation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *p, const G4Material *, G4double cut) override
 
virtual void PrintInfo () override
 
- 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 &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
void PrintInfoDefinition (const G4ParticleDefinition &part)
 
virtual void StartTracking (G4Track *) override
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double SampleSubCutSecondaries (std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual 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 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=nullptr)
 
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=nullptr)
 
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 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, G4bool lock=true)
 
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 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
 
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 *) override
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
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 ()
 

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 G4MuIonisation.hh.

Constructor & Destructor Documentation

G4MuIonisation::G4MuIonisation ( const G4String name = "muIoni")
explicit

Definition at line 99 of file G4MuIonisation.cc.

100  : G4VEnergyLossProcess(name),
101  theParticle(nullptr),
102  theBaseParticle(nullptr),
103  isInitialised(false)
104 {
105  mass = ratio = 0;
108 }
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetSecondaryParticle(const G4ParticleDefinition *p)
static G4Electron * Electron()
Definition: G4Electron.cc:94

Here is the call graph for this function:

G4MuIonisation::~G4MuIonisation ( )
virtual

Definition at line 112 of file G4MuIonisation.cc.

113 {}

Member Function Documentation

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

Implements G4VEnergyLossProcess.

Definition at line 135 of file G4MuIonisation.cc.

137 {
138  if(!isInitialised) {
139 
140  theParticle = part;
141  theBaseParticle = bpart;
142 
143  mass = theParticle->GetPDGMass();
144  G4double q = theParticle->GetPDGCharge();
145 
147  G4double elow = 0.2*MeV;
148  G4double emax = param->MaxKinEnergy();
149  G4double ehigh = std::min(1*GeV, emax);
150 
151  // Bragg peak model
152  if (!EmModel(1)) {
153  if(q > 0.0) { SetEmModel(new G4BraggModel(),1); }
154  else { SetEmModel(new G4ICRU73QOModel(),1); }
155  }
156  EmModel(1)->SetLowEnergyLimit(param->MinKinEnergy());
157  EmModel(1)->SetHighEnergyLimit(elow);
158  AddEmModel(1, EmModel(1), new G4IonFluctuations());
159 
160  // high energy fluctuation model
162 
163  // moderate energy model
164  if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
165  EmModel(2)->SetLowEnergyLimit(elow);
166  EmModel(2)->SetHighEnergyLimit(ehigh);
167  AddEmModel(2, EmModel(2), FluctModel());
168 
169  // high energy model
170  if(ehigh < emax) {
171  if (!EmModel(3)) { SetEmModel(new G4MuBetheBlochModel(),3); }
172  EmModel(3)->SetLowEnergyLimit(ehigh);
173  EmModel(3)->SetHighEnergyLimit(emax);
174  AddEmModel(3, EmModel(3), FluctModel());
175  }
176  ratio = electron_mass_c2/mass;
177  isInitialised = true;
178  }
179 }
G4double MaxKinEnergy() const
void SetFluctModel(G4VEmFluctuationModel *)
static constexpr double electron_mass_c2
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4VEmFluctuationModel * FluctModel()
G4double MinKinEnergy() const
static const G4double emax
G4double GetPDGMass() const
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
G4double GetPDGCharge() const
G4VEmModel * EmModel(G4int index=1) const

Here is the call graph for this function:

G4bool G4MuIonisation::IsApplicable ( const G4ParticleDefinition p)
overridevirtual

Implements G4VEnergyLossProcess.

Definition at line 117 of file G4MuIonisation.cc.

118 {
119  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV);
120 }
G4double GetPDGMass() const
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetPDGCharge() const

Here is the call graph for this function:

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 124 of file G4MuIonisation.cc.

127 {
128  G4double x = 0.5*cut/electron_mass_c2;
129  G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
130  return mass*(gam - 1.0);
131 }
static constexpr double electron_mass_c2
double G4double
Definition: G4Types.hh:76
void G4MuIonisation::PrintInfo ( )
overridevirtual

Implements G4VEnergyLossProcess.

Definition at line 183 of file G4MuIonisation.cc.

184 {}

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