Geant4 9.6.0
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=0, const G4String &nam="DNATransformElectronModel")
 
virtual ~G4DNATransformElectronModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetVerbose (int)
 
void SetEpsilonEnergy (G4double)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 51 of file G4DNATransformElectronModel.hh.

Constructor & Destructor Documentation

◆ G4DNATransformElectronModel()

G4DNATransformElectronModel::G4DNATransformElectronModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNATransformElectronModel" 
)

Definition at line 36 of file G4DNATransformElectronModel.cc.

37 :
38 G4VEmModel(nam),fIsInitialised(false)
39{
40 fVerboseLevel = 0 ;
41 SetLowEnergyLimit(0.*eV);
42 SetHighEnergyLimit(0.025*eV);
44 // fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45 fpWaterDensity = 0;
46 fpWaterDensity = 0;
47 fEpsilon = 0.0001*eV;
48}
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4DNATransformElectronModel()

G4DNATransformElectronModel::~G4DNATransformElectronModel ( )
virtual

Definition at line 51 of file G4DNATransformElectronModel.cc.

52{}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 83 of file G4DNATransformElectronModel.cc.

88{
89#if G4VERBOSE
90 if (fVerboseLevel > 1)
91 G4cout << "Calling CrossSectionPerVolume() of G4DNATransformElectronModel" << G4endl;
92#endif
93
94 if(ekin - fEpsilon > HighEnergyLimit())
95 {
96 return 0.0;
97 }
98
99 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
100
101 if(waterDensity!= 0.0)
102 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
103 {
104 if (ekin - fEpsilon <= HighEnergyLimit())
105 {
106 return DBL_MAX;
107 }
108 }
109
110 return 0.0 ;
111}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
size_t GetIndex() const
Definition: G4Material.hh:261
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
#define DBL_MAX
Definition: templates.hh:83

◆ Initialise()

void G4DNATransformElectronModel::Initialise ( const G4ParticleDefinition particleDefinition,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 55 of file G4DNATransformElectronModel.cc.

57{
58#ifdef G4VERBOSE
59 if (fVerboseLevel)
60 G4cout << "Calling G4DNATransformElectronModel::Initialise()" << G4endl;
61#endif
62
63 if (particleDefinition != G4Electron::ElectronDefinition())
64 {
65 G4ExceptionDescription exceptionDescription ;
66 exceptionDescription << "Attempting to calculate cross section for wrong particle";
67 G4Exception("G4DNATransformElectronModel::CrossSectionPerVolume","G4DNATransformElectronModel001",
68 FatalErrorInArgument,exceptionDescription);
69 return;
70 }
71
72 // Initialize water density pointer
74
75 if(!fIsInitialised)
76 {
77 fIsInitialised = true;
79 }
80}
@ FatalErrorInArgument
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 114 of file G4DNATransformElectronModel.cc.

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

◆ SetEpsilonEnergy()

void G4DNATransformElectronModel::SetEpsilonEnergy ( G4double  eps)
inline

Definition at line 99 of file G4DNATransformElectronModel.hh.

100{
101 fEpsilon = eps ;
102}

◆ SetVerbose()

void G4DNATransformElectronModel::SetVerbose ( int  flag)
inline

Definition at line 94 of file G4DNATransformElectronModel.hh.

95{
96 fVerboseLevel = flag ;
97}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNATransformElectronModel::fParticleChangeForGamma
protected

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