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

#include <G4DNATransformElectronModel.hh>

+ Inheritance diagram for G4DNATransformElectronModel:

Public Member Functions

 G4DNATransformElectronModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNATransformElectronModel")
 
 ~G4DNATransformElectronModel () override
 
G4DNATransformElectronModeloperator= (const G4DNATransformElectronModel &right)=delete
 
 G4DNATransformElectronModel (const G4DNATransformElectronModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetVerbose (int)
 
void SetEpsilonEnergy (G4double)
 
- 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 ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
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
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

When an electron reaches the highest energy domain of G4DNATransformElectronModel, it is then automatically converted into a solvated electron without thermalization displacement (assumed to be already thermalized).

Definition at line 59 of file G4DNATransformElectronModel.hh.

Constructor & Destructor Documentation

◆ G4DNATransformElectronModel() [1/2]

G4DNATransformElectronModel::G4DNATransformElectronModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNATransformElectronModel" )

Definition at line 36 of file G4DNATransformElectronModel.cc.

38 :
39 G4VEmModel(nam)
40{
41 fVerboseLevel = 0;
42 SetLowEnergyLimit(0. * eV);
43 SetHighEnergyLimit(0.025 * eV);
45 fpWaterDensity = nullptr;
46 fpWaterDensity = nullptr;
47 fEpsilon = 0.0001 * eV;
48}
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNATransformElectronModel()

G4DNATransformElectronModel::~G4DNATransformElectronModel ( )
overridedefault

◆ G4DNATransformElectronModel() [2/2]

G4DNATransformElectronModel::G4DNATransformElectronModel ( const G4DNATransformElectronModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNATransformElectronModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 86 of file G4DNATransformElectronModel.cc.

91{
92#if MODEL_VERBOSE
93 if (fVerboseLevel > 1)
94 G4cout << "Calling CrossSectionPerVolume() of G4DNATransformElectronModel" << G4endl;
95#endif
96
97 if(ekin - fEpsilon > HighEnergyLimit())
98 {
99 return 0.0;
100 }
101
102 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
103
104 if(waterDensity != 0.0)
105 {
106 return DBL_MAX;
107 }
108
109 return 0.0;
110}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t GetIndex() const
G4double HighEnergyLimit() const
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4DNATransformElectronModel::Initialise ( const G4ParticleDefinition * particleDefinition,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 55 of file G4DNATransformElectronModel.cc.

58{
59#ifdef MODEL_VERBOSE
60 if (fVerboseLevel)
61 G4cout << "Calling G4DNATransformElectronModel::Initialise()" << G4endl;
62#endif
63
64 if(particleDefinition->GetParticleName() != "e-")
65 {
67 errMsg << "Attempting to calculate cross section for wrong particle";
68 G4Exception("G4DNATransformElectronModel::CrossSectionPerVolume",
69 "G4DNATransformElectronModel001", FatalErrorInArgument, errMsg);
70 return;
71 }
72
73 // Initialize water density pointer
74 fpWaterDensity = G4DNAMolecularMaterial::Instance()->
75 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
76
77 if(!fIsInitialised)
78 {
79 fIsInitialised = true;
81 }
82}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ operator=()

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

◆ SampleSecondaries()

void G4DNATransformElectronModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * particle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 114 of file G4DNATransformElectronModel.cc.

119{
120#if MODEL_VERBOSE
121 if (fVerboseLevel)
122 G4cout << "Calling SampleSecondaries() of G4DNATransformElectronModel" << G4endl;
123#endif
124
125 G4double k = particle->GetKineticEnergy();
126
127// if (k - fEpsilon <= HighEnergyLimit()) // should be already checked
128// {
133// }
134}
@ fStopAndKill
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *pFinalPosition=nullptr)
static G4DNAChemistryManager * Instance()
G4double GetKineticEnergy() const
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetEpsilonEnergy()

void G4DNATransformElectronModel::SetEpsilonEnergy ( G4double eps)
inline

Definition at line 107 of file G4DNATransformElectronModel.hh.

108{
109 fEpsilon = eps ;
110}

◆ SetVerbose()

void G4DNATransformElectronModel::SetVerbose ( int flag)
inline

Definition at line 102 of file G4DNATransformElectronModel.hh.

103{
104 fVerboseLevel = flag ;
105}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNATransformElectronModel::fParticleChangeForGamma
protected

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