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

#include <G4BraggNoDeltaModel.hh>

+ Inheritance diagram for G4BraggNoDeltaModel:

Public Member Functions

 G4BraggNoDeltaModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggNoD")
 
 ~G4BraggNoDeltaModel () override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
G4BraggNoDeltaModeloperator= (const G4BraggNoDeltaModel &right)=delete
 
 G4BraggNoDeltaModel (const G4BraggNoDeltaModel &)=delete
 
- Public Member Functions inherited from G4BraggIonModel
 G4BraggIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
 
 ~G4BraggIonModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
 
void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss) override
 
G4BraggIonModeloperator= (const G4BraggIonModel &right)=delete
 
 G4BraggIonModel (const G4BraggIonModel &)=delete
 
- Public Member Functions inherited from G4BraggModel
 G4BraggModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")
 
 ~G4BraggModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
 
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
G4BraggModeloperator= (const G4BraggModel &right)=delete
 
 G4BraggModel (const G4BraggModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual void StartTracking (G4Track *)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4BraggModel
void SetParticle (const G4ParticleDefinition *p)
 
G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
void SetChargeSquareRatio (G4double val)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4BraggModel
const G4ParticleDefinitionparticle = nullptr
 
G4ParticleDefinitiontheElectron = nullptr
 
G4ParticleChangeForLossfParticleChange = nullptr
 
const G4MaterialcurrentMaterial = nullptr
 
const G4MaterialbaseMaterial = nullptr
 
G4EmCorrectionscorr = nullptr
 
G4double mass = 0.0
 
G4double spin = 0.0
 
G4double chargeSquare = 1.0
 
G4double massRate = 1.0
 
G4double ratio = 1.0
 
G4double protonMassAMU = 1.007276
 
G4double lowestKinEnergy
 
G4double theZieglerFactor
 
G4double expStopPower125
 
G4int iMolecula = -1
 
G4int iPSTAR = -1
 
G4int iICRU90 = -1
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 
- Static Protected Attributes inherited from G4BraggModel
static G4ICRU90StoppingDatafICRU90 = nullptr
 
static G4PSTARStoppingfPSTAR = nullptr
 

Detailed Description

Definition at line 53 of file G4BraggNoDeltaModel.hh.

Constructor & Destructor Documentation

◆ G4BraggNoDeltaModel() [1/2]

G4BraggNoDeltaModel::G4BraggNoDeltaModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "BraggNoD" )
explicit

Definition at line 53 of file G4BraggNoDeltaModel.cc.

54 :
55 G4BraggIonModel(p, nam)
56{}
G4BraggIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")

◆ ~G4BraggNoDeltaModel()

G4BraggNoDeltaModel::~G4BraggNoDeltaModel ( )
override

Definition at line 60 of file G4BraggNoDeltaModel.cc.

61{}

◆ G4BraggNoDeltaModel() [2/2]

G4BraggNoDeltaModel::G4BraggNoDeltaModel ( const G4BraggNoDeltaModel & )
delete

Member Function Documentation

◆ ComputeDEDXPerVolume()

G4double G4BraggNoDeltaModel::ComputeDEDXPerVolume ( const G4Material * material,
const G4ParticleDefinition * pd,
G4double kineticEnergy,
G4double cutEnergy )
overridevirtual

Reimplemented from G4BraggIonModel.

Definition at line 65 of file G4BraggNoDeltaModel.cc.

69{
70 return
71 G4BraggIonModel::ComputeDEDXPerVolume(material, pd, kinEnergy, DBL_MAX);
72}
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
#define DBL_MAX
Definition templates.hh:62

◆ CrossSectionPerVolume()

G4double G4BraggNoDeltaModel::CrossSectionPerVolume ( const G4Material * ,
const G4ParticleDefinition * ,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4BraggIonModel.

Definition at line 76 of file G4BraggNoDeltaModel.cc.

80{
81 return 0.0;
82}

◆ operator=()

G4BraggNoDeltaModel & G4BraggNoDeltaModel::operator= ( const G4BraggNoDeltaModel & right)
delete

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