Geant4 9.6.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)
 
G4double PartialCrossSection (G4double energy, G4int level)
 
- 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

Definition at line 42 of file G4DNAEmfietzoglouExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAEmfietzoglouExcitationModel()

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

Definition at line 41 of file G4DNAEmfietzoglouExcitationModel.cc.

43 :G4VEmModel(nam),isInitialised(false)
44{
45 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
46 fpWaterDensity = 0;
47
48 lowEnergyLimit = 8.23 * eV;
49 highEnergyLimit = 10 * MeV;
50 SetLowEnergyLimit(lowEnergyLimit);
51 SetHighEnergyLimit(highEnergyLimit);
52
53 nLevels = waterExcitation.NumberOfLevels();
54
55 verboseLevel= 0;
56 // Verbosity scale:
57 // 0 = nothing
58 // 1 = warning for energy non-conservation
59 // 2 = details of energy budget
60 // 3 = calculation of cross sections, file openings, sampling of atoms
61 // 4 = entering in methods
62
63 if (verboseLevel > 3)
64 if( verboseLevel>0 )
65 {
66 G4cout << "Emfietzoglou Excitation model is constructed " << G4endl
67 << "Energy range: "
68 << lowEnergyLimit / eV << " eV - "
69 << highEnergyLimit / MeV << " MeV"
70 << G4endl;
71 }
73}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4DNAEmfietzoglouExcitationModel()

G4DNAEmfietzoglouExcitationModel::~G4DNAEmfietzoglouExcitationModel ( )
virtual

Definition at line 77 of file G4DNAEmfietzoglouExcitationModel.cc.

78{}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 126 of file G4DNAEmfietzoglouExcitationModel.cc.

131{
132 if (verboseLevel > 3)
133 G4cout << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouExcitationModel" << G4endl;
134
135 // Calculate total cross section for model
136
137 G4double sigma=0;
138
139 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
140
141 if(waterDensity!= 0.0)
142 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
143 {
144
145 if (particleDefinition == G4Electron::ElectronDefinition())
146 {
147 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
148 {
149 sigma = Sum(ekin);
150 }
151 }
152
153 if (verboseLevel > 2)
154 {
155 G4cout << "__________________________________" << G4endl;
156 G4cout << "°°° G4DNAEmfietzoglouExcitationModel - XS INFO START" << G4endl;
157 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
158 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
159 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
160 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
161 G4cout << "°°° G4DNAEmfietzoglouExcitationModel - XS INFO END" << G4endl;
162 }
163
164 }
165
166 return sigma*material->GetAtomicNumDensityVector()[1];
167}
double G4double
Definition: G4Types.hh:64
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
size_t GetIndex() const
Definition: G4Material.hh:261

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 82 of file G4DNAEmfietzoglouExcitationModel.cc.

84{
85
86 if (verboseLevel > 3)
87 G4cout << "Calling G4DNAEmfietzoglouExcitationModel::Initialise()" << G4endl;
88
89 // Energy limits
90
91 if (LowEnergyLimit() < lowEnergyLimit)
92 {
93 G4cout << "G4DNAEmfietzoglouExcitationModel: low energy limit increased from " <<
94 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
95 SetLowEnergyLimit(lowEnergyLimit);
96 }
97
98 if (HighEnergyLimit() > highEnergyLimit)
99 {
100 G4cout << "G4DNAEmfietzoglouExcitationModel: high energy limit decreased from " <<
101 HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
102 SetHighEnergyLimit(highEnergyLimit);
103 }
104
105 //
106 if( verboseLevel>0 )
107 {
108 G4cout << "Emfietzoglou Excitation model is initialized " << G4endl
109 << "Energy range: "
110 << LowEnergyLimit() / eV << " eV - "
111 << HighEnergyLimit() / MeV << " MeV"
112 << G4endl;
113 }
114
115 // Initialize water density pointer
117
118 if (isInitialised) { return; }
120 isInitialised = true;
121
122}
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ PartialCrossSection()

G4double G4DNAEmfietzoglouExcitationModel::PartialCrossSection ( G4double  energy,
G4int  level 
)

Definition at line 203 of file G4DNAEmfietzoglouExcitationModel.cc.

204{
205 // Aj T
206 // Sigma(T) = ------------- (Bj / T) ln(Cj ---) [1 - Bj / T]^Pj
207 // 2 pi alpha0 R
208 //
209 // Sigma is the macroscopic cross section = N sigma, where N = number of target particles per unit volume
210 // and sigma is the microscopic cross section
211 // T is the incoming electron kinetic energy
212 // alpha0 is the Bohr Radius (Bohr_radius)
213 // Aj, Bj, Cj & Pj are parameters that can be found in Emfietzoglou's papers
214 //
215 // From Phys. Med. Biol. 48 (2003) 2355-2371, D.Emfietzoglou,
216 // Monte Carlo Simulation of the energy loss of low energy electrons in liquid Water
217 //
218 // Scaling for macroscopic cross section: number of water moleculs per unit volume
219 // const G4double sigma0 = (10. / 3.343e22) * cm2;
220
221 const G4double density = 3.34192e+19 * mm3;
222
223 const G4double aj[]={0.0205, 0.0209, 0.0130, 0.0026, 0.0025};
224 const G4double cj[]={4.9801, 3.3850, 2.8095, 1.9242, 3.4624};
225 const G4double pj[]={0.4757, 0.3483, 0.4443, 0.3429, 0.4379};
226 const G4double r = 13.6 * eV;
227
228 G4double sigma = 0.;
229
230 G4double exc = waterExcitation.ExcitationEnergy(level);
231
232 if (t >= exc)
233 {
234 G4double excitationSigma = ( aj[level] / (2.*pi*Bohr_radius))
235 * (exc / t)
236 * std::log(cj[level]*(t/r))
237 * std::pow((1.- (exc/t)), pj[level]);
238 sigma = excitationSigma / density;
239 }
240
241 return sigma;
242}
const G4double pi

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 171 of file G4DNAEmfietzoglouExcitationModel.cc.

176{
177
178 if (verboseLevel > 3)
179 G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouExcitationModel" << G4endl;
180
181 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
182
183 G4int level = RandomSelect(electronEnergy0);
184
185 G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
186 G4double newEnergy = electronEnergy0 - excitationEnergy;
187
188 if (electronEnergy0 < highEnergyLimit)
189 {
193
194 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
196 level,
197 theIncomingTrack);
198 }
199}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:66
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAEmfietzoglouExcitationModel::fParticleChangeForGamma
protected

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