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

#include <G4DNAEmfietzoglouExcitationModel.hh>

+ Inheritance diagram for G4DNAEmfietzoglouExcitationModel:

Public Member Functions

 G4DNAEmfietzoglouExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouExcitationModel")
 
virtual ~G4DNAEmfietzoglouExcitationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new 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 SelectStationary (G4bool input)
 
- 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 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

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

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

Definition at line 50 of file G4DNAEmfietzoglouExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAEmfietzoglouExcitationModel()

G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAEmfietzoglouExcitationModel" 
)

Definition at line 47 of file G4DNAEmfietzoglouExcitationModel.cc.

49:G4VEmModel(nam),isInitialised(false)
50{
51 fpMolWaterDensity = 0;
52
53 verboseLevel= 0;
54 // Verbosity scale:
55 // 0 = nothing
56 // 1 = warning for energy non-conservation
57 // 2 = details of energy budget
58 // 3 = calculation of cross sections, file openings, sampling of atoms
59 // 4 = entering in methods
60
61 if( verboseLevel>0 )
62 {
63 G4cout << "Emfietzoglou excitation model is constructed " << G4endl;
64 }
66
67 SetLowEnergyLimit(8.*eV);
68 SetHighEnergyLimit(10.*keV);
69
70 // Selection of stationary mode
71 statCode = false;
72}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ ~G4DNAEmfietzoglouExcitationModel()

G4DNAEmfietzoglouExcitationModel::~G4DNAEmfietzoglouExcitationModel ( )
virtual

Definition at line 76 of file G4DNAEmfietzoglouExcitationModel.cc.

77{
78 // Cross section
79
80 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
81 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
82 {
83 G4DNACrossSectionDataSet* table = pos->second;
84 delete table;
85 }
86
87}

Member Function Documentation

◆ CrossSectionPerVolume()

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

sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;

Reimplemented from G4VEmModel.

Definition at line 141 of file G4DNAEmfietzoglouExcitationModel.cc.

146{
147 if (verboseLevel > 3)
148 G4cout << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouExcitationModel" << G4endl;
149
150 if (particleDefinition != G4Electron::ElectronDefinition()) return 0;
151
152 // Calculate total cross section for model
153
154 G4double sigma=0;
155
156 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
157
158 const G4String& particleName = particleDefinition->GetParticleName();
159
160 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
161 {
162 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
163 pos = tableData.find(particleName);
164
165 if (pos != tableData.end())
166 {
167 G4DNACrossSectionDataSet* table = pos->second;
168 if (table != 0) sigma = table->FindValue(ekin);
169 }
170 else
171 {
172 G4Exception("G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume","em0002",
173 FatalException,"Model not applicable to particle type.");
174 }
175 }
176
177 if (verboseLevel > 2)
178 {
179 G4cout << "__________________________________" << G4endl;
180 G4cout << "G4DNAEmfietzoglouExcitationModel - XS INFO START" << G4endl;
181 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
182 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
183 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
184 //G4cout << " Cross section per water molecule (cm^-1)=" <<
185 ///sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
186 G4cout << "G4DNAEmfietzoglouExcitationModel - XS INFO END" << G4endl;
187 }
188
189 return sigma*waterDensity;
190}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
size_t GetIndex() const
Definition: G4Material.hh:258
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ Initialise()

void G4DNAEmfietzoglouExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 91 of file G4DNAEmfietzoglouExcitationModel.cc.

93{
94
95 if (verboseLevel > 3)
96 G4cout << "Calling G4DNAEmfietzoglouExcitationModel::Initialise()" << G4endl;
97
98 G4String fileElectron("dna/sigma_excitation_e_emfietzoglou");
99
101
103
104 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
105
106 // *** ELECTRON
107
108 electron = electronDef->GetParticleName();
109
110 tableFile[electron] = fileElectron;
111
112 // Cross section
113
115 tableE->LoadData(fileElectron);
116
117 tableData[electron] = tableE;
118
119 //
120
121 if( verboseLevel>0 )
122 {
123 G4cout << "Emfietzoglou excitation model is initialized " << G4endl
124 << "Energy range: "
125 << LowEnergyLimit() / eV << " eV - "
126 << HighEnergyLimit() / keV << " keV for "
127 << particle->GetParticleName()
128 << G4endl;
129 }
130
131 // Initialize water density pointer
133
134 if (isInitialised) return;
136 isInitialised = true;
137}
virtual G4bool LoadData(const G4String &argFileName)
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 194 of file G4DNAEmfietzoglouExcitationModel.cc.

199{
200
201 if (verboseLevel > 3)
202 G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouExcitationModel" << G4endl;
203
204 G4double k = aDynamicParticle->GetKineticEnergy();
205
206 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
207
208 G4int level = RandomSelect(k,particleName);
209 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
210 G4double newEnergy = k - excitationEnergy;
211
212 if (newEnergy > 0)
213 {
215
216 if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
218
220 }
221
222 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
224 level,
225 theIncomingTrack);
226}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNAEmfietzoglouExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 120 of file G4DNAEmfietzoglouExcitationModel.hh.

121{
122 statCode = input;
123}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAEmfietzoglouExcitationModel::fParticleChangeForGamma
protected

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