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

#include <G4hImpactIonisation.hh>

+ Inheritance diagram for G4hImpactIonisation:

Public Member Functions

 G4hImpactIonisation (const G4String &processName="hImpactIoni")
 
 ~G4hImpactIonisation ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, enum G4ForceCondition *condition)
 
void PrintInfoDefinition () const
 
void SetHighEnergyForProtonParametrisation (G4double energy)
 
void SetLowEnergyForProtonParametrisation (G4double energy)
 
void SetHighEnergyForAntiProtonParametrisation (G4double energy)
 
void SetLowEnergyForAntiProtonParametrisation (G4double energy)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
void SetElectronicStoppingPowerModel (const G4ParticleDefinition *aParticle, const G4String &dedxTable)
 
void SetNuclearStoppingPowerModel (const G4String &dedxTable)
 
void SetNuclearStoppingOn ()
 
void SetNuclearStoppingOff ()
 
void SetBarkasOn ()
 
void SetBarkasOff ()
 
void SetPixe (const G4bool)
 
G4VParticleChangeAlongStepDoIt (const G4Track &trackData, const G4Step &stepData)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)
 
G4double ComputeDEDX (const G4ParticleDefinition *aParticle, const G4MaterialCutsCouple *couple, G4double kineticEnergy)
 
void SetCutForSecondaryPhotons (G4double cut)
 
void SetCutForAugerElectrons (G4double cut)
 
void ActivateAugerElectronProduction (G4bool val)
 
void SetPixeCrossSectionK (const G4String &name)
 
void SetPixeCrossSectionL (const G4String &name)
 
void SetPixeCrossSectionM (const G4String &name)
 
void SetPixeProjectileMinEnergy (G4double energy)
 
void SetPixeProjectileMaxEnergy (G4double energy)
 
- Public Member Functions inherited from G4hRDEnergyLoss
 G4hRDEnergyLoss (const G4String &)
 
 ~G4hRDEnergyLoss ()
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, enum G4ForceCondition *condition)=0
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)=0
 
- 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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4hRDEnergyLoss
static G4int GetNumberOfProcesses ()
 
static void SetNumberOfProcesses (G4int number)
 
static void PlusNumberOfProcesses ()
 
static void MinusNumberOfProcesses ()
 
static void SetdRoverRange (G4double value)
 
static void SetRndmStep (G4bool value)
 
static void SetEnlossFluc (G4bool value)
 
static void SetStepFunction (G4double c1, G4double c2)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4hRDEnergyLoss
G4bool CutsWhereModified ()
 
- 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 ()
 
- Static Protected Member Functions inherited from G4hRDEnergyLoss
static void BuildDEDXTable (const G4ParticleDefinition &aParticleType)
 
- Protected Attributes inherited from G4hRDEnergyLoss
const G4double MaxExcitationNumber
 
const G4double probLimFluct
 
const long nmaxDirectFluct
 
const long nmaxCont1
 
const long nmaxCont2
 
G4PhysicsTabletheLossTable
 
G4double linLossLimit
 
G4double MinKineticEnergy
 
- 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
 
- Static Protected Attributes inherited from G4hRDEnergyLoss
static G4ThreadLocal G4PhysicsTabletheDEDXpTable = 0
 
static G4ThreadLocal G4PhysicsTabletheDEDXpbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheRangepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheRangepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheInverseRangepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheInverseRangepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheLabTimepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheLabTimepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheProperTimepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheProperTimepbarTable = 0
 
static G4ThreadLocal G4PhysicsTable ** RecorderOfpProcess = 0
 
static G4ThreadLocal G4PhysicsTable ** RecorderOfpbarProcess = 0
 
static G4ThreadLocal G4int CounterOfpProcess = 0
 
static G4ThreadLocal G4int CounterOfpbarProcess = 0
 
static G4ThreadLocal G4double ParticleMass
 
static G4ThreadLocal G4double ptableElectronCutInRange = 0.0
 
static G4ThreadLocal G4double pbartableElectronCutInRange = 0.0
 
static G4ThreadLocal G4double Charge
 
static G4ThreadLocal G4double LowestKineticEnergy = 1e-05
 
static G4ThreadLocal G4double HighestKineticEnergy = 1.e5
 
static G4ThreadLocal G4int TotBin = 360
 
static G4ThreadLocal G4double RTable
 
static G4ThreadLocal G4double LOGRTable
 
static G4ThreadLocal G4double dRoverRange = 0.20
 
static G4ThreadLocal G4double finalRange = 0.2
 
static G4ThreadLocal G4double c1lim = 0.20
 
static G4ThreadLocal G4double c2lim = 0.32
 
static G4ThreadLocal G4double c3lim = -0.032
 
static G4ThreadLocal G4bool rndmStepFlag = false
 
static G4ThreadLocal G4bool EnlossFlucFlag = true
 

Detailed Description

Definition at line 74 of file G4hImpactIonisation.hh.

Constructor & Destructor Documentation

◆ G4hImpactIonisation()

G4hImpactIonisation::G4hImpactIonisation ( const G4String processName = "hImpactIoni")

Definition at line 82 of file G4hImpactIonisation.cc.

83 : G4hRDEnergyLoss(processName),
84 betheBlochModel(0),
85 protonModel(0),
86 antiprotonModel(0),
87 theIonEffChargeModel(0),
88 theNuclearStoppingModel(0),
89 theIonChuFluctuationModel(0),
90 theIonYangFluctuationModel(0),
91 protonTable("ICRU_R49p"),
92 antiprotonTable("ICRU_R49p"),
93 theNuclearTable("ICRU_R49"),
94 nStopping(true),
95 theBarkas(true),
96 theMeanFreePathTable(0),
97 paramStepLimit (0.005),
98 pixeCrossSectionHandler(0)
99{
100 InitializeMe();
101}

◆ ~G4hImpactIonisation()

G4hImpactIonisation::~G4hImpactIonisation ( )

Definition at line 131 of file G4hImpactIonisation.cc.

132{
133 if (theMeanFreePathTable)
134 {
135 theMeanFreePathTable->clearAndDestroy();
136 delete theMeanFreePathTable;
137 }
138
139 if (betheBlochModel) delete betheBlochModel;
140 if (protonModel) delete protonModel;
141 if (antiprotonModel) delete antiprotonModel;
142 if (theNuclearStoppingModel) delete theNuclearStoppingModel;
143 if (theIonEffChargeModel) delete theIonEffChargeModel;
144 if (theIonChuFluctuationModel) delete theIonChuFluctuationModel;
145 if (theIonYangFluctuationModel) delete theIonYangFluctuationModel;
146
147 delete pixeCrossSectionHandler;
148
149 // ---- MGP ---- The following is to be checked
150 // if (shellVacancy) delete shellVacancy;
151
152 cutForDelta.clear();
153}
void clearAndDestroy()

Member Function Documentation

◆ ActivateAugerElectronProduction()

void G4hImpactIonisation::ActivateAugerElectronProduction ( G4bool  val)

Definition at line 1706 of file G4hImpactIonisation.cc.

1707{
1708 atomicDeexcitation.ActivateAugerElectronProduction(val);
1709}
void ActivateAugerElectronProduction(G4bool val)
Set threshold energy for Auger electron production.

◆ AlongStepDoIt()

G4VParticleChange * G4hImpactIonisation::AlongStepDoIt ( const G4Track trackData,
const G4Step stepData 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 717 of file G4hImpactIonisation.cc.

719{
720 // compute the energy loss after a step
723 G4double finalT = 0.;
724
726
727 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
728 const G4Material* material = couple->GetMaterial();
729
730 // get the actual (true) Step length from step
731 const G4double stepLength = step.GetStepLength() ;
732
733 const G4DynamicParticle* particle = track.GetDynamicParticle() ;
734
735 G4double kineticEnergy = particle->GetKineticEnergy() ;
736 G4double massRatio = proton_mass_c2/(particle->GetMass()) ;
737 G4double tScaled = kineticEnergy * massRatio ;
738 G4double eLoss = 0.0 ;
739 G4double nLoss = 0.0 ;
740
741
742 // very small particle energy
743 if (kineticEnergy < MinKineticEnergy)
744 {
745 eLoss = kineticEnergy ;
746 // particle energy outside tabulated energy range
747 }
748
749 else if( kineticEnergy > HighestKineticEnergy)
750 {
751 eLoss = stepLength * fdEdx ;
752 // big step
753 }
754 else if (stepLength >= fRangeNow )
755 {
756 eLoss = kineticEnergy ;
757
758 // tabulated range
759 }
760 else
761 {
762 // step longer than linear step limit
763 if(stepLength > linLossLimit * fRangeNow)
764 {
765 G4double rScaled = fRangeNow * massRatio * chargeSquare ;
766 G4double sScaled = stepLength * massRatio * chargeSquare ;
767
768 if(charge > 0.0)
769 {
770 eLoss = G4EnergyLossTables::GetPreciseEnergyFromRange(proton,rScaled, couple) -
771 G4EnergyLossTables::GetPreciseEnergyFromRange(proton,rScaled-sScaled,couple) ;
772
773 }
774 else
775 {
776 // Antiproton
777 eLoss = G4EnergyLossTables::GetPreciseEnergyFromRange(antiproton,rScaled,couple) -
778 G4EnergyLossTables::GetPreciseEnergyFromRange(antiproton,rScaled-sScaled,couple) ;
779 }
780 eLoss /= massRatio ;
781
782 // Barkas correction at big step
783 eLoss += fBarkas * stepLength;
784
785 // step shorter than linear step limit
786 }
787 else
788 {
789 eLoss = stepLength *fdEdx ;
790 }
791 if (nStopping && tScaled < protonHighEnergy)
792 {
793 nLoss = (theNuclearStoppingModel->TheValue(particle, material)) * stepLength;
794 }
795 }
796
797 if (eLoss < 0.0) eLoss = 0.0;
798
799 finalT = kineticEnergy - eLoss - nLoss;
800
801 if ( EnlossFlucFlag && 0.0 < eLoss && finalT > MinKineticEnergy)
802 {
803
804 // now the electron loss with fluctuation
805 eLoss = ElectronicLossFluctuation(particle, couple, eLoss, stepLength) ;
806 if (eLoss < 0.0) eLoss = 0.0;
807 finalT = kineticEnergy - eLoss - nLoss;
808 }
809
810 // stop particle if the kinetic energy <= MinKineticEnergy
811 if (finalT*massRatio <= MinKineticEnergy )
812 {
813
814 finalT = 0.0;
817 else
819 }
820
822 eLoss = kineticEnergy-finalT;
823
825 return &aParticleChange ;
826}
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4double GetPreciseEnergyFromRange(const G4ParticleDefinition *aParticle, G4double range, const G4Material *aMaterial)
const G4Material * GetMaterial() const
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
G4ProcessManager * GetProcessManager() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual G4double TheValue(const G4DynamicParticle *particle, const G4Material *material)=0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4double MinKineticEnergy
static G4ThreadLocal G4double HighestKineticEnergy
static G4ThreadLocal G4bool EnlossFlucFlag

◆ BuildPhysicsTable()

void G4hImpactIonisation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 190 of file G4hImpactIonisation.cc.

193{
194
195 // Verbose print-out
196 if(verboseLevel > 0)
197 {
198 G4cout << "G4hImpactIonisation::BuildPhysicsTable for "
199 << particleDef.GetParticleName()
200 << " mass(MeV)= " << particleDef.GetPDGMass()/MeV
201 << " charge= " << particleDef.GetPDGCharge()/eplus
202 << " type= " << particleDef.GetParticleType()
203 << G4endl;
204
205 if(verboseLevel > 1)
206 {
207 G4ProcessVector* pv = particleDef.GetProcessManager()->GetProcessList();
208
209 G4cout << " 0: " << (*pv)[0]->GetProcessName() << " " << (*pv)[0]
210 << " 1: " << (*pv)[1]->GetProcessName() << " " << (*pv)[1]
211 // << " 2: " << (*pv)[2]->GetProcessName() << " " << (*pv)[2]
212 << G4endl;
213 G4cout << "ionModel= " << theIonEffChargeModel
214 << " MFPtable= " << theMeanFreePathTable
215 << " iniMass= " << initialMass
216 << G4endl;
217 }
218 }
219 // End of verbose print-out
220
221 if (particleDef.GetParticleType() == "nucleus" &&
222 particleDef.GetParticleName() != "GenericIon" &&
223 particleDef.GetParticleSubType() == "generic")
224 {
225
226 G4EnergyLossTables::Register(&particleDef,
233 proton_mass_c2/particleDef.GetPDGMass(),
234 TotBin);
235
236 return;
237 }
238
239 if( !CutsWhereModified() && theLossTable) return;
240
241 InitializeParametrisation() ;
244
245 charge = particleDef.GetPDGCharge() / eplus;
246 chargeSquare = charge*charge ;
247
248 const G4ProductionCutsTable* theCoupleTable=
250 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
251
252 cutForDelta.clear();
253 cutForGamma.clear();
254
255 for (G4int j=0; j<numOfCouples; ++j) {
256
257 // get material parameters needed for the energy loss calculation
258 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
259 const G4Material* material= couple->GetMaterial();
260
261 // the cut cannot be below lowest limit
262 G4double tCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[j];
264
265 G4double excEnergy = material->GetIonisation()->GetMeanExcitationEnergy();
266
267 tCut = std::max(tCut,excEnergy);
268 cutForDelta.push_back(tCut);
269
270 // the cut cannot be below lowest limit
271 tCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[j];
273 tCut = std::max(tCut,minGammaEnergy);
274 cutForGamma.push_back(tCut);
275 }
276
277 if(verboseLevel > 0) {
278 G4cout << "Cuts are defined " << G4endl;
279 }
280
281 if(0.0 < charge)
282 {
283 {
284 BuildLossTable(*proton) ;
285
286 // The following vector has a fixed dimension (see src/G4hImpactLoss.cc for more details)
287 // It happended in the past that caused memory corruption errors. The problem is still pending, even if temporary solved
288 // G4cout << "[NOTE]: __LINE__=" << __LINE__ << ", particleDef=" << particleDef.GetParticleName() << ", proton=" << proton << ", theLossTable=" << theLossTable << ", CounterOfpProcess=" << CounterOfpProcess << G4endl;
289
292 }
293 } else {
294 {
295 BuildLossTable(*antiproton) ;
296
297 // The following vector has a fixed dimension (see src/G4hImpactLoss.cc for more details)
298 // It happended in the past that caused memory corruption errors. The problem is still pending, even if temporary solved
299 // G4cout << "[NOTE]: __LINE__=" << __LINE__ << ", particleDef=" << particleDef.GetParticleName() << ", antiproton=" << antiproton << ", theLossTable=" << theLossTable << ", CounterOfpbarProcess=" << CounterOfpbarProcess << G4endl;
300
303 }
304 }
305
306 if(verboseLevel > 0) {
307 G4cout << "G4hImpactIonisation::BuildPhysicsTable: "
308 << "Loss table is built "
309 // << theLossTable
310 << G4endl;
311 }
312
313 BuildLambdaTable(particleDef) ;
314 // BuildDataForFluorescence(particleDef);
315
316 if(verboseLevel > 1) {
317 G4cout << (*theMeanFreePathTable) << G4endl;
318 }
319
320 if(verboseLevel > 0) {
321 G4cout << "G4hImpactIonisation::BuildPhysicsTable: "
322 << "DEDX table will be built "
323 // << theDEDXpTable << " " << theDEDXpbarTable
324 // << " " << theRangepTable << " " << theRangepbarTable
325 << G4endl;
326 }
327
328 BuildDEDXTable(particleDef) ;
329
330 if(verboseLevel > 1) {
331 G4cout << (*theDEDXpTable) << G4endl;
332 }
333
334 if((&particleDef == proton) || (&particleDef == antiproton)) PrintInfoDefinition() ;
335
336 if(verboseLevel > 0) {
337 G4cout << "G4hImpactIonisation::BuildPhysicsTable: end for "
338 << particleDef.GetParticleName() << G4endl;
339 }
340}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static void Register(const G4ParticleDefinition *p, const G4PhysicsTable *tDEDX, const G4PhysicsTable *tRange, const G4PhysicsTable *tInverseRange, const G4PhysicsTable *tLabTime, const G4PhysicsTable *tProperTime, G4double lowestKineticEnergy, G4double highestKineticEnergy, G4double massRatio, G4int NumberOfBins)
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4int verboseLevel
Definition: G4VProcess.hh:360
static G4ThreadLocal G4PhysicsTable * theDEDXpTable
static G4ThreadLocal G4PhysicsTable * theInverseRangepTable
static G4ThreadLocal G4PhysicsTable ** RecorderOfpProcess
static G4ThreadLocal G4int TotBin
static G4ThreadLocal G4PhysicsTable * theProperTimepTable
G4PhysicsTable * theLossTable
static G4ThreadLocal G4PhysicsTable * theLabTimepTable
static G4ThreadLocal G4PhysicsTable * theRangepTable
static G4ThreadLocal G4int CounterOfpProcess
static G4ThreadLocal G4int CounterOfpbarProcess
static G4ThreadLocal G4PhysicsTable ** RecorderOfpbarProcess
static G4ThreadLocal G4double LowestKineticEnergy
G4bool CutsWhereModified()
static void BuildDEDXTable(const G4ParticleDefinition &aParticleType)

◆ ComputeDEDX()

G4double G4hImpactIonisation::ComputeDEDX ( const G4ParticleDefinition aParticle,
const G4MaterialCutsCouple couple,
G4double  kineticEnergy 
)

Definition at line 1265 of file G4hImpactIonisation.cc.

1268{
1269 const G4Material* material = couple->GetMaterial();
1271 G4AntiProton* antiproton = G4AntiProton::AntiProton();
1272 G4double dedx = 0.;
1273
1274 G4double tScaled = kineticEnergy * proton_mass_c2 / (aParticle->GetPDGMass()) ;
1275 charge = aParticle->GetPDGCharge() ;
1276
1277 if( charge > 0.)
1278 {
1279 if (tScaled > protonHighEnergy)
1280 {
1281 dedx = G4EnergyLossTables::GetDEDX(proton,tScaled,couple) ;
1282 }
1283 else
1284 {
1285 dedx = ProtonParametrisedDEDX(couple,tScaled) ;
1286 }
1287 }
1288 else
1289 {
1290 if (tScaled > antiprotonHighEnergy)
1291 {
1292 dedx = G4EnergyLossTables::GetDEDX(antiproton,tScaled,couple);
1293 }
1294 else
1295 {
1296 dedx = AntiProtonParametrisedDEDX(couple,tScaled) ;
1297 }
1298 }
1299 dedx *= theIonEffChargeModel->TheValue(aParticle, material, kineticEnergy) ;
1300
1301 return dedx ;
1302}
static G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
G4double GetPDGCharge() const

◆ GetContinuousStepLimit()

G4double G4hImpactIonisation::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
inlinevirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 293 of file G4hImpactIonisation.hh.

297{
298 G4double step = GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
299
300 // ---- MGP ---- The following line, taken as is from G4hLowEnergyIonisation,
301 // is meaningless: currentMinimumStep is passed by value,
302 // therefore any local modification to it has no effect
303
304 if ((step > 0.) && (step < currentMinimumStep)) currentMinimumStep = step ;
305
306 return step ;
307}
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

◆ GetMeanFreePath()

G4double G4hImpactIonisation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
enum G4ForceCondition condition 
)
virtual

Implements G4hRDEnergyLoss.

Definition at line 588 of file G4hImpactIonisation.cc.

591{
592 const G4DynamicParticle* dynamicParticle = track.GetDynamicParticle();
593 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
594 const G4Material* material = couple->GetMaterial();
595
596 G4double meanFreePath = DBL_MAX;
597 // ---- MGP ---- What is the meaning of the local variable isOutOfRange?
598 G4bool isOutRange = false;
599
601
602 G4double kineticEnergy = (dynamicParticle->GetKineticEnergy())*initialMass/(dynamicParticle->GetMass());
603 charge = dynamicParticle->GetCharge()/eplus;
604 chargeSquare = theIonEffChargeModel->TheValue(dynamicParticle, material);
605
606 if (kineticEnergy < LowestKineticEnergy)
607 {
608 meanFreePath = DBL_MAX;
609 }
610 else
611 {
612 if (kineticEnergy > HighestKineticEnergy) kineticEnergy = HighestKineticEnergy;
613 meanFreePath = (((*theMeanFreePathTable)(couple->GetIndex()))->
614 GetValue(kineticEnergy,isOutRange))/chargeSquare;
615 }
616
617 return meanFreePath ;
618}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
bool G4bool
Definition: G4Types.hh:86
G4double GetCharge() const
#define DBL_MAX
Definition: templates.hh:62

◆ IsApplicable()

G4bool G4hImpactIonisation::IsApplicable ( const G4ParticleDefinition particle)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 310 of file G4hImpactIonisation.hh.

311{
312 // ---- MGP ---- Better criterion for applicability to be defined;
313 // now hard-coded particle mass > 0.1 * proton_mass
314
315 return (particle.GetPDGCharge() != 0.0 && particle.GetPDGMass() > CLHEP::proton_mass_c2*0.1);
316}

◆ PostStepDoIt()

G4VParticleChange * G4hImpactIonisation::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
virtual

Implements G4hRDEnergyLoss.

Definition at line 951 of file G4hImpactIonisation.cc.

953{
954 // Units are expressed in GEANT4 internal units.
955
956 // std::cout << "----- Calling PostStepDoIt ----- " << std::endl;
957
959 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
960 const G4DynamicParticle* aParticle = track.GetDynamicParticle() ;
961
962 // Some kinematics
963
964 G4ParticleDefinition* definition = track.GetDefinition();
965 G4double mass = definition->GetPDGMass();
966 G4double kineticEnergy = aParticle->GetKineticEnergy();
967 G4double totalEnergy = kineticEnergy + mass ;
968 G4double pSquare = kineticEnergy *( totalEnergy + mass) ;
969 G4double eSquare = totalEnergy * totalEnergy;
970 G4double betaSquare = pSquare / eSquare;
971 G4ThreeVector particleDirection = aParticle->GetMomentumDirection() ;
972
973 G4double gamma = kineticEnergy / mass + 1.;
974 G4double r = electron_mass_c2 / mass;
975 G4double tMax = 2. * electron_mass_c2 *(gamma*gamma - 1.) / (1. + 2.*gamma*r + r*r);
976
977 // Validity range for delta electron cross section
978 G4double deltaCut = cutForDelta[couple->GetIndex()];
979
980 // This should not be a case
981 if (deltaCut >= tMax)
983
984 G4double xc = deltaCut / tMax;
985 G4double rate = tMax / totalEnergy;
986 rate = rate*rate ;
987 G4double spin = aParticle->GetDefinition()->GetPDGSpin() ;
988
989 // Sampling follows ...
990 G4double x = 0.;
991 G4double gRej = 0.;
992
993 do {
994 x = xc / (1. - (1. - xc) * G4UniformRand());
995
996 if (0.0 == spin)
997 {
998 gRej = 1.0 - betaSquare * x ;
999 }
1000 else if (0.5 == spin)
1001 {
1002 gRej = (1. - betaSquare * x + 0.5 * x*x * rate) / (1. + 0.5 * rate) ;
1003 }
1004 else
1005 {
1006 gRej = (1. - betaSquare * x ) * (1. + x/(3.*xc)) +
1007 x*x * rate * (1. + 0.5*x/xc) / 3.0 /
1008 (1. + 1./(3.*xc) + rate *(1.+ 0.5/xc) / 3.);
1009 }
1010
1011 } while( G4UniformRand() > gRej );
1012
1013 G4double deltaKineticEnergy = x * tMax;
1014 G4double deltaTotalMomentum = std::sqrt(deltaKineticEnergy *
1015 (deltaKineticEnergy + 2. * electron_mass_c2 ));
1016 G4double totalMomentum = std::sqrt(pSquare) ;
1017 G4double cosTheta = deltaKineticEnergy * (totalEnergy + electron_mass_c2) / (deltaTotalMomentum*totalMomentum);
1018
1019 // protection against cosTheta > 1 or < -1
1020 if ( cosTheta < -1. ) cosTheta = -1.;
1021 if ( cosTheta > 1. ) cosTheta = 1.;
1022
1023 // direction of the delta electron
1024 G4double phi = twopi * G4UniformRand();
1025 G4double sinTheta = std::sqrt(1. - cosTheta*cosTheta);
1026 G4double dirX = sinTheta * std::cos(phi);
1027 G4double dirY = sinTheta * std::sin(phi);
1028 G4double dirZ = cosTheta;
1029
1030 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
1031 deltaDirection.rotateUz(particleDirection);
1032
1033 // create G4DynamicParticle object for delta ray
1034 G4DynamicParticle* deltaRay = new G4DynamicParticle;
1035 deltaRay->SetKineticEnergy( deltaKineticEnergy );
1036 deltaRay->SetMomentumDirection(deltaDirection.x(),
1037 deltaDirection.y(),
1038 deltaDirection.z());
1040
1041 // fill aParticleChange
1042 G4double finalKineticEnergy = kineticEnergy - deltaKineticEnergy;
1043 std::size_t totalNumber = 1;
1044
1045 // Atomic relaxation
1046
1047 // ---- MGP ---- Temporary limitation: currently PIXE only for incident protons
1048
1049 std::size_t nSecondaries = 0;
1050 std::vector<G4DynamicParticle*>* secondaryVector = 0;
1051
1052 if (definition == G4Proton::ProtonDefinition())
1053 {
1054 const G4Material* material = couple->GetMaterial();
1055
1056 // Lazy initialization of pixeCrossSectionHandler
1057 if (pixeCrossSectionHandler == 0)
1058 {
1059 // Instantiate pixeCrossSectionHandler with selected shell cross section models
1060 // Ownership of interpolation is transferred to pixeCrossSectionHandler
1061 G4IInterpolator* interpolation = new G4LogLogInterpolator();
1062 pixeCrossSectionHandler =
1063 new G4PixeCrossSectionHandler(interpolation,modelK,modelL,modelM,eMinPixe,eMaxPixe);
1064 G4String fileName("proton");
1065 pixeCrossSectionHandler->LoadShellData(fileName);
1066 // pixeCrossSectionHandler->PrintData();
1067 }
1068
1069 // Select an atom in the current material based on the total shell cross sections
1070 G4int Z = pixeCrossSectionHandler->SelectRandomAtom(material,kineticEnergy);
1071 // std::cout << "G4hImpactIonisation::PostStepDoIt - Z = " << Z << std::endl;
1072
1073 // G4double microscopicCross = MicroscopicCrossSection(*definition,
1074 // kineticEnergy,
1075 // Z, deltaCut);
1076 // G4double crossFromShells = pixeCrossSectionHandler->FindValue(Z,kineticEnergy);
1077
1078 //std::cout << "G4hImpactIonisation: Z= "
1079 // << Z
1080 // << ", energy = "
1081 // << kineticEnergy/MeV
1082 // <<" MeV, microscopic = "
1083 // << microscopicCross/barn
1084 // << " barn, from shells = "
1085 // << crossFromShells/barn
1086 // << " barn"
1087 // << std::endl;
1088
1089 // Select a shell in the target atom based on the individual shell cross sections
1090 G4int shellIndex = pixeCrossSectionHandler->SelectRandomShell(Z,kineticEnergy);
1091
1093 const G4AtomicShell* atomicShell = transitionManager->Shell(Z,shellIndex);
1094 G4double bindingEnergy = atomicShell->BindingEnergy();
1095
1096 // if (verboseLevel > 1)
1097 // {
1098 // G4cout << "G4hImpactIonisation::PostStepDoIt - Z = "
1099 // << Z
1100 // << ", shell = "
1101 // << shellIndex
1102 // << ", bindingE (keV) = "
1103 // << bindingEnergy/keV
1104 // << G4endl;
1105 // }
1106
1107 // Generate PIXE if binding energy larger than cut for photons or electrons
1108
1109 G4ParticleDefinition* type = 0;
1110
1111 if (finalKineticEnergy >= bindingEnergy)
1112 // && (bindingEnergy >= minGammaEnergy || bindingEnergy >= minElectronEnergy) )
1113 {
1114 // Vacancy in subshell shellIndex; shellId is the subshell identifier in EADL jargon
1115 G4int shellId = atomicShell->ShellId();
1116 // Atomic relaxation: generate secondaries
1117 secondaryVector = atomicDeexcitation.GenerateParticles(Z, shellId);
1118
1119 // ---- Debug ----
1120 //std::cout << "ShellId = "
1121 // <<shellId << " ---- Atomic relaxation secondaries: ---- "
1122 // << secondaryVector->size()
1123 // << std::endl;
1124
1125 // ---- End debug ---
1126
1127 if (secondaryVector != 0)
1128 {
1129 nSecondaries = secondaryVector->size();
1130 for (std::size_t i = 0; i<nSecondaries; i++)
1131 {
1132 G4DynamicParticle* aSecondary = (*secondaryVector)[i];
1133 if (aSecondary)
1134 {
1135 G4double e = aSecondary->GetKineticEnergy();
1136 type = aSecondary->GetDefinition();
1137
1138 // ---- Debug ----
1139 //if (type == G4Gamma::GammaDefinition())
1140 // {
1141 // std::cout << "Z = " << Z
1142 // << ", shell: " << shellId
1143 // << ", PIXE photon energy (keV) = " << e/keV
1144 // << std::endl;
1145 // }
1146 // ---- End debug ---
1147
1148 if (e < finalKineticEnergy &&
1149 ((type == G4Gamma::Gamma() && e > minGammaEnergy ) ||
1150 (type == G4Electron::Electron() && e > minElectronEnergy )))
1151 {
1152 // Subtract the energy of the emitted secondary from the primary
1153 finalKineticEnergy -= e;
1154 totalNumber++;
1155 // ---- Debug ----
1156 //if (type == G4Gamma::GammaDefinition())
1157 // {
1158 // std::cout << "Z = " << Z
1159 // << ", shell: " << shellId
1160 // << ", PIXE photon energy (keV) = " << e/keV
1161 // << std::endl;
1162 // }
1163 // ---- End debug ---
1164 }
1165 else
1166 {
1167 // The atomic relaxation product has energy below the cut
1168 // ---- Debug ----
1169 // if (type == G4Gamma::GammaDefinition())
1170 // {
1171 // std::cout << "Z = " << Z
1172 //
1173 // << ", PIXE photon energy = " << e/keV
1174 // << " keV below threshold " << minGammaEnergy/keV << " keV"
1175 // << std::endl;
1176 // }
1177 // ---- End debug ---
1178
1179 delete aSecondary;
1180 (*secondaryVector)[i] = 0;
1181 }
1182 }
1183 }
1184 }
1185 }
1186 }
1187
1188
1189 // Save delta-electrons
1190
1191 G4double eDeposit = 0.;
1192
1193 if (finalKineticEnergy > MinKineticEnergy)
1194 {
1195 G4double finalPx = totalMomentum*particleDirection.x() - deltaTotalMomentum*deltaDirection.x();
1196 G4double finalPy = totalMomentum*particleDirection.y() - deltaTotalMomentum*deltaDirection.y();
1197 G4double finalPz = totalMomentum*particleDirection.z() - deltaTotalMomentum*deltaDirection.z();
1198 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz) ;
1199 finalPx /= finalMomentum;
1200 finalPy /= finalMomentum;
1201 finalPz /= finalMomentum;
1202
1203 aParticleChange.ProposeMomentumDirection( finalPx,finalPy,finalPz );
1204 }
1205 else
1206 {
1207 eDeposit = finalKineticEnergy;
1208 finalKineticEnergy = 0.;
1209 aParticleChange.ProposeMomentumDirection(particleDirection.x(),
1210 particleDirection.y(),
1211 particleDirection.z());
1212 if(!aParticle->GetDefinition()->GetProcessManager()->
1213 GetAtRestProcessVector()->size())
1215 else
1217 }
1218
1219 aParticleChange.ProposeEnergy(finalKineticEnergy);
1222 aParticleChange.AddSecondary(deltaRay);
1223
1224 // ---- Debug ----
1225 // std::cout << "RDHadronIonisation - finalKineticEnergy (MeV) = "
1226 // << finalKineticEnergy/MeV
1227 // << ", delta KineticEnergy (keV) = "
1228 // << deltaKineticEnergy/keV
1229 // << ", energy deposit (MeV) = "
1230 // << eDeposit/MeV
1231 // << std::endl;
1232 // ---- End debug ---
1233
1234 // Save Fluorescence and Auger
1235
1236 if (secondaryVector != 0)
1237 {
1238 for (std::size_t l = 0; l < nSecondaries; l++)
1239 {
1240 G4DynamicParticle* secondary = (*secondaryVector)[l];
1241 if (secondary) aParticleChange.AddSecondary(secondary);
1242
1243 // ---- Debug ----
1244 //if (secondary != 0)
1245 // {
1246 // if (secondary->GetDefinition() == G4Gamma::GammaDefinition())
1247 // {
1248 // G4double eX = secondary->GetKineticEnergy();
1249 // std::cout << " PIXE photon of energy " << eX/keV
1250 // << " keV added to ParticleChange; total number of secondaries is " << totalNumber
1251 // << std::endl;
1252 // }
1253 //}
1254 // ---- End debug ---
1255
1256 }
1257 delete secondaryVector;
1258 }
1259
1261}
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
G4double BindingEnergy() const
G4int ShellId() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
static G4AtomicTransitionManager * Instance()
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void AddSecondary(G4Track *aSecondary)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4int SelectRandomShell(G4int Z, G4double e) const
void LoadShellData(const G4String &dataFile)
G4int SelectRandomAtom(const G4Material *material, G4double e) const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
G4ParticleDefinition * GetDefinition() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double bindingEnergy(G4int A, G4int Z)

◆ PrintInfoDefinition()

void G4hImpactIonisation::PrintInfoDefinition ( ) const

Definition at line 1713 of file G4hImpactIonisation.cc.

1714{
1715 G4String comments = " Knock-on electron cross sections . ";
1716 comments += "\n Good description above the mean excitation energy.\n";
1717 comments += " Delta ray energy sampled from differential Xsection.";
1718
1719 G4cout << G4endl << GetProcessName() << ": " << comments
1720 << "\n PhysicsTables from " << LowestKineticEnergy / eV << " eV "
1721 << " to " << HighestKineticEnergy / TeV << " TeV "
1722 << " in " << TotBin << " bins."
1723 << "\n Electronic stopping power model is "
1724 << protonTable
1725 << "\n from " << protonLowEnergy / keV << " keV "
1726 << " to " << protonHighEnergy / MeV << " MeV " << "." << G4endl ;
1727 G4cout << "\n Parametrisation model for antiprotons is "
1728 << antiprotonTable
1729 << "\n from " << antiprotonLowEnergy / keV << " keV "
1730 << " to " << antiprotonHighEnergy / MeV << " MeV " << "." << G4endl ;
1731 if(theBarkas){
1732 G4cout << " Parametrization of the Barkas effect is switched on."
1733 << G4endl ;
1734 }
1735 if(nStopping) {
1736 G4cout << " Nuclear stopping power model is " << theNuclearTable
1737 << G4endl ;
1738 }
1739
1740 G4bool printHead = true;
1741
1742 const G4ProductionCutsTable* theCoupleTable=
1744 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
1745
1746 // loop for materials
1747
1748 for (G4int j=0 ; j < numOfCouples; ++j) {
1749
1750 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
1751 const G4Material* material= couple->GetMaterial();
1752 G4double deltaCutNow = cutForDelta[(couple->GetIndex())] ;
1753 G4double excitationEnergy = material->GetIonisation()->GetMeanExcitationEnergy();
1754
1755 if(excitationEnergy > deltaCutNow) {
1756 if(printHead) {
1757 printHead = false ;
1758
1759 G4cout << " material min.delta energy(keV) " << G4endl;
1760 G4cout << G4endl;
1761 }
1762
1763 G4cout << std::setw(20) << material->GetName()
1764 << std::setw(15) << excitationEnergy/keV << G4endl;
1765 }
1766 }
1767}
const G4String & GetName() const
Definition: G4Material.hh:172
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

Referenced by BuildPhysicsTable().

◆ SetBarkasOff()

void G4hImpactIonisation::SetBarkasOff ( )
inline

Definition at line 147 of file G4hImpactIonisation.hh.

147{theBarkas = false;};

◆ SetBarkasOn()

void G4hImpactIonisation::SetBarkasOn ( )
inline

Definition at line 144 of file G4hImpactIonisation.hh.

144{theBarkas = true;};

◆ SetCutForAugerElectrons()

void G4hImpactIonisation::SetCutForAugerElectrons ( G4double  cut)

Definition at line 1699 of file G4hImpactIonisation.cc.

1700{
1701 minElectronEnergy = cut;
1702}

◆ SetCutForSecondaryPhotons()

void G4hImpactIonisation::SetCutForSecondaryPhotons ( G4double  cut)

Definition at line 1692 of file G4hImpactIonisation.cc.

1693{
1694 minGammaEnergy = cut;
1695}

◆ SetElectronicStoppingPowerModel()

void G4hImpactIonisation::SetElectronicStoppingPowerModel ( const G4ParticleDefinition aParticle,
const G4String dedxTable 
)

Definition at line 156 of file G4hImpactIonisation.cc.

159{
160 if (particle->GetPDGCharge() > 0 )
161 {
162 // Positive charge
163 SetProtonElectronicStoppingPowerModel(dedxTable) ;
164 }
165 else
166 {
167 // Antiprotons
168 SetAntiProtonElectronicStoppingPowerModel(dedxTable) ;
169 }
170}

◆ SetHighEnergyForAntiProtonParametrisation()

void G4hImpactIonisation::SetHighEnergyForAntiProtonParametrisation ( G4double  energy)
inline

Definition at line 109 of file G4hImpactIonisation.hh.

109{antiprotonHighEnergy = energy;} ;
G4double energy(const ThreeVector &p, const G4double m)

◆ SetHighEnergyForProtonParametrisation()

void G4hImpactIonisation::SetHighEnergyForProtonParametrisation ( G4double  energy)
inline

Definition at line 99 of file G4hImpactIonisation.hh.

99{protonHighEnergy = energy;} ;

◆ SetLowEnergyForAntiProtonParametrisation()

void G4hImpactIonisation::SetLowEnergyForAntiProtonParametrisation ( G4double  energy)
inline

Definition at line 114 of file G4hImpactIonisation.hh.

114{antiprotonLowEnergy = energy;} ;

◆ SetLowEnergyForProtonParametrisation()

void G4hImpactIonisation::SetLowEnergyForProtonParametrisation ( G4double  energy)
inline

Definition at line 104 of file G4hImpactIonisation.hh.

104{protonLowEnergy = energy;} ;

◆ SetNuclearStoppingOff()

void G4hImpactIonisation::SetNuclearStoppingOff ( )
inline

Definition at line 141 of file G4hImpactIonisation.hh.

141{nStopping = false;};

◆ SetNuclearStoppingOn()

void G4hImpactIonisation::SetNuclearStoppingOn ( )
inline

Definition at line 138 of file G4hImpactIonisation.hh.

138{nStopping = true;};

Referenced by SetNuclearStoppingPowerModel().

◆ SetNuclearStoppingPowerModel()

void G4hImpactIonisation::SetNuclearStoppingPowerModel ( const G4String dedxTable)
inline

Definition at line 130 of file G4hImpactIonisation.hh.

131 {theNuclearTable = dedxTable; SetNuclearStoppingOn();};

◆ SetPixe()

void G4hImpactIonisation::SetPixe ( const  G4bool)
inline

Definition at line 150 of file G4hImpactIonisation.hh.

150{pixeIsActive = true;};

◆ SetPixeCrossSectionK()

void G4hImpactIonisation::SetPixeCrossSectionK ( const G4String name)
inline

Definition at line 176 of file G4hImpactIonisation.hh.

176{ modelK = name; }
const char * name(G4int ptype)

◆ SetPixeCrossSectionL()

void G4hImpactIonisation::SetPixeCrossSectionL ( const G4String name)
inline

Definition at line 177 of file G4hImpactIonisation.hh.

177{ modelL = name; }

◆ SetPixeCrossSectionM()

void G4hImpactIonisation::SetPixeCrossSectionM ( const G4String name)
inline

Definition at line 178 of file G4hImpactIonisation.hh.

178{ modelM = name; }

◆ SetPixeProjectileMaxEnergy()

void G4hImpactIonisation::SetPixeProjectileMaxEnergy ( G4double  energy)
inline

Definition at line 180 of file G4hImpactIonisation.hh.

180{ eMaxPixe = energy; }

◆ SetPixeProjectileMinEnergy()

void G4hImpactIonisation::SetPixeProjectileMinEnergy ( G4double  energy)
inline

Definition at line 179 of file G4hImpactIonisation.hh.

179{ eMinPixe = energy; }

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