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

#include <G4DNATransformElectronModel.hh>

Inheritance diagram for G4DNATransformElectronModel:
Collaboration 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 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

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

Constructor & Destructor Documentation

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

Definition at line 38 of file G4DNATransformElectronModel.cc.

39  :
40  G4VEmModel(nam),fIsInitialised(false)
41 {
42  fVerboseLevel = 0;
43  SetLowEnergyLimit(0. * eV);
44  SetHighEnergyLimit(0.025 * eV);
46  fpWaterDensity = 0;
47  fpWaterDensity = 0;
48  fEpsilon = 0.0001 * eV;
49 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * fParticleChangeForGamma
static constexpr double eV
Definition: G4SIunits.hh:215
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739

Here is the call graph for this function:

G4DNATransformElectronModel::~G4DNATransformElectronModel ( )
virtual

Definition at line 52 of file G4DNATransformElectronModel.cc.

53 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 87 of file G4DNATransformElectronModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 57 of file G4DNATransformElectronModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 115 of file G4DNATransformElectronModel.cc.

120 {
121 #if MODEL_VERBOSE
122  if (fVerboseLevel)
123  G4cout << "Calling SampleSecondaries() of G4DNATransformElectronModel" << G4endl;
124 #endif
125 
126  G4double k = particle->GetKineticEnergy();
127 
128 // if (k - fEpsilon <= HighEnergyLimit()) // should be already checked
129 // {
134 // }
135 }
G4double GetKineticEnergy() const
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAChemistryManager * Instance()
const G4Track * GetCurrentTrack() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

void G4DNATransformElectronModel::SetEpsilonEnergy ( G4double  eps)
inline

Definition at line 108 of file G4DNATransformElectronModel.hh.

109 {
110  fEpsilon = eps ;
111 }
static const G4double eps
void G4DNATransformElectronModel::SetVerbose ( int  flag)
inline

Definition at line 103 of file G4DNATransformElectronModel.hh.

104 {
105  fVerboseLevel = flag ;
106 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNATransformElectronModel::fParticleChangeForGamma
protected

Definition at line 87 of file G4DNATransformElectronModel.hh.


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