Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4hImpactIonisation Class Reference

#include <G4hImpactIonisation.hh>

Inheritance diagram for G4hImpactIonisation:
Collaboration 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 ()
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
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)
 
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 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
 
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
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
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
 
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
 
- Static Protected Attributes inherited from G4hRDEnergyLoss
static G4ThreadLocal
G4PhysicsTable
theDEDXpTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theDEDXpbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theRangepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theRangepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theInverseRangepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theInverseRangepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theLabTimepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theLabTimepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theProperTimepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theProperTimepbarTable = 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 75 of file G4hImpactIonisation.hh.

Constructor & Destructor Documentation

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

Definition at line 83 of file G4hImpactIonisation.cc.

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

Definition at line 132 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

Member Function Documentation

void G4hImpactIonisation::ActivateAugerElectronProduction ( G4bool  val)

Definition at line 1707 of file G4hImpactIonisation.cc.

1708 {
1709  atomicDeexcitation.ActivateAugerElectronProduction(val);
1710 }
void ActivateAugerElectronProduction(G4bool val)

Here is the call graph for this function:

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 718 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

void G4hImpactIonisation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 191 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

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

Definition at line 1266 of file G4hImpactIonisation.cc.

1269 {
1270  const G4Material* material = couple->GetMaterial();
1272  G4AntiProton* antiproton = G4AntiProton::AntiProton();
1273  G4double dedx = 0.;
1274 
1275  G4double tScaled = kineticEnergy * proton_mass_c2 / (aParticle->GetPDGMass()) ;
1276  charge = aParticle->GetPDGCharge() ;
1277 
1278  if( charge > 0.)
1279  {
1280  if (tScaled > protonHighEnergy)
1281  {
1282  dedx = G4EnergyLossTables::GetDEDX(proton,tScaled,couple) ;
1283  }
1284  else
1285  {
1286  dedx = ProtonParametrisedDEDX(couple,tScaled) ;
1287  }
1288  }
1289  else
1290  {
1291  if (tScaled > antiprotonHighEnergy)
1292  {
1293  dedx = G4EnergyLossTables::GetDEDX(antiproton,tScaled,couple);
1294  }
1295  else
1296  {
1297  dedx = AntiProtonParametrisedDEDX(couple,tScaled) ;
1298  }
1299  }
1300  dedx *= theIonEffChargeModel->TheValue(aParticle, material, kineticEnergy) ;
1301 
1302  return dedx ;
1303 }
string material
Definition: eplot.py:19
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
G4double GetPDGMass() const
virtual G4double TheValue(const G4DynamicParticle *particle, const G4Material *material)=0
static G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Implements G4VContinuousDiscreteProcess.

Definition at line 294 of file G4hImpactIonisation.hh.

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

Here is the call graph for this function:

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

Implements G4hRDEnergyLoss.

Definition at line 589 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

G4bool G4hImpactIonisation::IsApplicable ( const G4ParticleDefinition particle)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 311 of file G4hImpactIonisation.hh.

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

Here is the call graph for this function:

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

Implements G4hRDEnergyLoss.

Definition at line 952 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

void G4hImpactIonisation::PrintInfoDefinition ( ) const

Definition at line 1714 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4hImpactIonisation::SetBarkasOff ( )
inline

Definition at line 148 of file G4hImpactIonisation.hh.

148 {theBarkas = false;};
void G4hImpactIonisation::SetBarkasOn ( )
inline

Definition at line 145 of file G4hImpactIonisation.hh.

145 {theBarkas = true;};
void G4hImpactIonisation::SetCutForAugerElectrons ( G4double  cut)

Definition at line 1700 of file G4hImpactIonisation.cc.

1701 {
1702  minElectronEnergy = cut;
1703 }

Here is the caller graph for this function:

void G4hImpactIonisation::SetCutForSecondaryPhotons ( G4double  cut)

Definition at line 1693 of file G4hImpactIonisation.cc.

1694 {
1695  minGammaEnergy = cut;
1696 }

Here is the caller graph for this function:

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

Definition at line 157 of file G4hImpactIonisation.cc.

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

Here is the call graph for this function:

void G4hImpactIonisation::SetHighEnergyForAntiProtonParametrisation ( G4double  energy)
inline

Definition at line 110 of file G4hImpactIonisation.hh.

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

Here is the call graph for this function:

void G4hImpactIonisation::SetHighEnergyForProtonParametrisation ( G4double  energy)
inline

Definition at line 100 of file G4hImpactIonisation.hh.

100 {protonHighEnergy = energy;} ;
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void G4hImpactIonisation::SetLowEnergyForAntiProtonParametrisation ( G4double  energy)
inline

Definition at line 115 of file G4hImpactIonisation.hh.

115 {antiprotonLowEnergy = energy;} ;
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void G4hImpactIonisation::SetLowEnergyForProtonParametrisation ( G4double  energy)
inline

Definition at line 105 of file G4hImpactIonisation.hh.

105 {protonLowEnergy = energy;} ;
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void G4hImpactIonisation::SetNuclearStoppingOff ( )
inline

Definition at line 142 of file G4hImpactIonisation.hh.

142 {nStopping = false;};
void G4hImpactIonisation::SetNuclearStoppingOn ( )
inline

Definition at line 139 of file G4hImpactIonisation.hh.

139 {nStopping = true;};

Here is the caller graph for this function:

void G4hImpactIonisation::SetNuclearStoppingPowerModel ( const G4String dedxTable)
inline

Definition at line 131 of file G4hImpactIonisation.hh.

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

Here is the call graph for this function:

void G4hImpactIonisation::SetPixe ( const G4bool  )
inline

Definition at line 151 of file G4hImpactIonisation.hh.

151 {pixeIsActive = true;};
void G4hImpactIonisation::SetPixeCrossSectionK ( const G4String name)
inline

Definition at line 177 of file G4hImpactIonisation.hh.

177 { modelK = name; }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

void G4hImpactIonisation::SetPixeCrossSectionL ( const G4String name)
inline

Definition at line 178 of file G4hImpactIonisation.hh.

178 { modelL = name; }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

void G4hImpactIonisation::SetPixeCrossSectionM ( const G4String name)
inline

Definition at line 179 of file G4hImpactIonisation.hh.

179 { modelM = name; }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

void G4hImpactIonisation::SetPixeProjectileMaxEnergy ( G4double  energy)
inline

Definition at line 181 of file G4hImpactIonisation.hh.

181 { eMaxPixe = energy; }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4hImpactIonisation::SetPixeProjectileMinEnergy ( G4double  energy)
inline

Definition at line 180 of file G4hImpactIonisation.hh.

180 { eMinPixe = energy; }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

Here is the caller graph for this function:


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