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

#include <G4DNAMillerGreenExcitationModel.hh>

+ Inheritance diagram for G4DNAMillerGreenExcitationModel:

Public Member Functions

 G4DNAMillerGreenExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAMillerGreenExcitationModel")
 
virtual ~G4DNAMillerGreenExcitationModel ()
 
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)
 
- 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 43 of file G4DNAMillerGreenExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMillerGreenExcitationModel()

G4DNAMillerGreenExcitationModel::G4DNAMillerGreenExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAMillerGreenExcitationModel" 
)

Definition at line 41 of file G4DNAMillerGreenExcitationModel.cc.

43 :G4VEmModel(nam),isInitialised(false)
44{
45 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
46 fpMolWaterDensity = 0;
47
48 nLevels=0;
49 kineticEnergyCorrection[0]=0.;
50 kineticEnergyCorrection[1]=0.;
51 kineticEnergyCorrection[2]=0.;
52 kineticEnergyCorrection[3]=0.;
53
54 verboseLevel= 0;
55 // Verbosity scale:
56 // 0 = nothing
57 // 1 = warning for energy non-conservation
58 // 2 = details of energy budget
59 // 3 = calculation of cross sections, file openings, sampling of atoms
60 // 4 = entering in methods
61
62 if( verboseLevel>0 )
63 {
64 G4cout << "Miller & Green excitation model is constructed " << G4endl;
65 }
67}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma

◆ ~G4DNAMillerGreenExcitationModel()

G4DNAMillerGreenExcitationModel::~G4DNAMillerGreenExcitationModel ( )
virtual

Definition at line 71 of file G4DNAMillerGreenExcitationModel.cc.

72{}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 222 of file G4DNAMillerGreenExcitationModel.cc.

227{
228 if (verboseLevel > 3)
229 G4cout << "Calling CrossSectionPerVolume() of G4DNAMillerGreenExcitationModel" << G4endl;
230
231 // Calculate total cross section for model
232
233 G4DNAGenericIonsManager *instance;
235
236 if (
237 particleDefinition != G4Proton::ProtonDefinition()
238 &&
239 particleDefinition != instance->GetIon("hydrogen")
240 &&
241 particleDefinition != instance->GetIon("alpha++")
242 &&
243 particleDefinition != instance->GetIon("alpha+")
244 &&
245 particleDefinition != instance->GetIon("helium")
246 )
247
248 return 0;
249
250 G4double lowLim = 0;
251 G4double highLim = 0;
252 G4double crossSection = 0.;
253
254 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
255
256 if(waterDensity!= 0.0)
257 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
258 {
259// G4cout << "WATER DENSITY = " << waterDensity*G4Material::GetMaterial("G4_WATER")->GetMassOfMolecule()/(g/cm3)
260// << G4endl;
261 const G4String& particleName = particleDefinition->GetParticleName();
262
263 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
264 pos1 = lowEnergyLimit.find(particleName);
265
266 if (pos1 != lowEnergyLimit.end())
267 {
268 lowLim = pos1->second;
269 }
270
271 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
272 pos2 = highEnergyLimit.find(particleName);
273
274 if (pos2 != highEnergyLimit.end())
275 {
276 highLim = pos2->second;
277 }
278
279 if (k >= lowLim && k < highLim)
280 {
281 crossSection = Sum(k,particleDefinition);
282
283 // add ONE or TWO electron-water excitation for alpha+ and helium
284 /*
285 if ( particleDefinition == instance->GetIon("alpha+")
286 ||
287 particleDefinition == instance->GetIon("helium")
288 )
289 {
290
291 G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
292 excitationXS->Initialise(G4Electron::ElectronDefinition());
293
294 G4double sigmaExcitation=0;
295 G4double tmp =0.;
296
297 if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation =
298 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
299 /material->GetAtomicNumDensityVector()[1];
300
301 if ( particleDefinition == instance->GetIon("alpha+") )
302 crossSection = crossSection + sigmaExcitation ;
303
304 if ( particleDefinition == instance->GetIon("helium") )
305 crossSection = crossSection + 2*sigmaExcitation ;
306
307 delete excitationXS;
308
309 // Alternative excitation model
310
311 G4DNABornExcitationModel * excitationXS = new G4DNABornExcitationModel();
312 excitationXS->Initialise(G4Electron::ElectronDefinition());
313
314 G4double sigmaExcitation=0;
315 G4double tmp=0;
316
317 if (k*0.511/3728 > 9*eV && k*0.511/3728 < 1*MeV ) sigmaExcitation =
318 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
319 /material->GetAtomicNumDensityVector()[1];
320
321 if ( particleDefinition == instance->GetIon("alpha+") )
322 crossSection = crossSection + sigmaExcitation ;
323
324 if ( particleDefinition == instance->GetIon("helium") )
325 crossSection = crossSection + 2*sigmaExcitation ;
326
327 delete excitationXS;
328
329 }
330*/
331
332 }
333
334 if (verboseLevel > 2)
335 {
336 G4cout << "__________________________________" << G4endl;
337 G4cout << "°°° G4DNAMillerGreenExcitationModel - XS INFO START" << G4endl;
338 G4cout << "°°° Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
339 G4cout << "°°° Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
340 G4cout << "°°° Cross section per water molecule (cm^-1)=" << crossSection*waterDensity/(1./cm) << G4endl;
341 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
342 G4cout << "°°° G4DNAMillerGreenExcitationModel - XS INFO END" << G4endl;
343 }
344 }
345 else
346 {
347 if (verboseLevel > 2)
348 {
349 G4cout << "MillerGreenExcitationModel : WARNING Water density is NULL" << G4endl;
350 }
351 }
352
353 return crossSection*waterDensity;
354 // return crossSection*material->GetAtomicNumDensityVector()[1];
355
356}
double G4double
Definition: G4Types.hh:64
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
size_t GetIndex() const
Definition: G4Material.hh:261
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88

◆ Initialise()

void G4DNAMillerGreenExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 76 of file G4DNAMillerGreenExcitationModel.cc.

78{
79
80 if (verboseLevel > 3)
81 G4cout << "Calling G4DNAMillerGreenExcitationModel::Initialise()" << G4endl;
82
83 // Energy limits
84
88 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
89 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
90 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
91 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
92
94 G4String hydrogen;
95 G4String alphaPlusPlus;
96 G4String alphaPlus;
97 G4String helium;
98
99 // LIMITS AND CONSTANTS
100
101 proton = protonDef->GetParticleName();
102 lowEnergyLimit[proton] = 10. * eV;
103 highEnergyLimit[proton] = 500. * keV;
104
105 kineticEnergyCorrection[0] = 1.;
106 slaterEffectiveCharge[0][0] = 0.;
107 slaterEffectiveCharge[1][0] = 0.;
108 slaterEffectiveCharge[2][0] = 0.;
109 sCoefficient[0][0] = 0.;
110 sCoefficient[1][0] = 0.;
111 sCoefficient[2][0] = 0.;
112
113 hydrogen = hydrogenDef->GetParticleName();
114 lowEnergyLimit[hydrogen] = 10. * eV;
115 highEnergyLimit[hydrogen] = 500. * keV;
116
117 kineticEnergyCorrection[0] = 1.;
118 slaterEffectiveCharge[0][0] = 0.;
119 slaterEffectiveCharge[1][0] = 0.;
120 slaterEffectiveCharge[2][0] = 0.;
121 sCoefficient[0][0] = 0.;
122 sCoefficient[1][0] = 0.;
123 sCoefficient[2][0] = 0.;
124
125 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
126 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
127 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
128
129 kineticEnergyCorrection[1] = 0.9382723/3.727417;
130 slaterEffectiveCharge[0][1]=0.;
131 slaterEffectiveCharge[1][1]=0.;
132 slaterEffectiveCharge[2][1]=0.;
133 sCoefficient[0][1]=0.;
134 sCoefficient[1][1]=0.;
135 sCoefficient[2][1]=0.;
136
137 alphaPlus = alphaPlusDef->GetParticleName();
138 lowEnergyLimit[alphaPlus] = 1. * keV;
139 highEnergyLimit[alphaPlus] = 400. * MeV;
140
141 kineticEnergyCorrection[2] = 0.9382723/3.727417;
142 slaterEffectiveCharge[0][2]=2.0;
143
144 // Following values provided by M. Dingfelder
145 slaterEffectiveCharge[1][2]=2.00;
146 slaterEffectiveCharge[2][2]=2.00;
147 //
148 sCoefficient[0][2]=0.7;
149 sCoefficient[1][2]=0.15;
150 sCoefficient[2][2]=0.15;
151
152 helium = heliumDef->GetParticleName();
153 lowEnergyLimit[helium] = 1. * keV;
154 highEnergyLimit[helium] = 400. * MeV;
155
156 kineticEnergyCorrection[3] = 0.9382723/3.727417;
157 slaterEffectiveCharge[0][3]=1.7;
158 slaterEffectiveCharge[1][3]=1.15;
159 slaterEffectiveCharge[2][3]=1.15;
160 sCoefficient[0][3]=0.5;
161 sCoefficient[1][3]=0.25;
162 sCoefficient[2][3]=0.25;
163
164 //
165
166 if (particle==protonDef)
167 {
168 SetLowEnergyLimit(lowEnergyLimit[proton]);
169 SetHighEnergyLimit(highEnergyLimit[proton]);
170 }
171
172 if (particle==hydrogenDef)
173 {
174 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
175 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
176 }
177
178 if (particle==alphaPlusPlusDef)
179 {
180 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
181 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
182 }
183
184 if (particle==alphaPlusDef)
185 {
186 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
187 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
188 }
189
190 if (particle==heliumDef)
191 {
192 SetLowEnergyLimit(lowEnergyLimit[helium]);
193 SetHighEnergyLimit(highEnergyLimit[helium]);
194 }
195
196 //
197
198 nLevels = waterExcitation.NumberOfLevels();
199
200 //
201 if( verboseLevel>0 )
202 {
203 G4cout << "Miller & Green excitation model is initialized " << G4endl
204 << "Energy range: "
205 << LowEnergyLimit() / eV << " eV - "
206 << HighEnergyLimit() / keV << " keV for "
207 << particle->GetParticleName()
208 << G4endl;
209 }
210
211 // Initialize water density pointer
213
214 if (isInitialised) { return; }
216 isInitialised = true;
217
218}
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 360 of file G4DNAMillerGreenExcitationModel.cc.

365{
366
367 if (verboseLevel > 3)
368 G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
369
370 G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
371
372 G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
373
374 // G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
375
376 // Dingfelder's excitation levels
377 const G4double excitation[]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
378 G4double excitationEnergy = excitation[level];
379
380 G4double newEnergy = particleEnergy0 - excitationEnergy;
381
382 if (newEnergy>0)
383 {
387
388 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
390 level,
391 theIncomingTrack);
392
393 }
394
395}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:66
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)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMillerGreenExcitationModel::fParticleChangeForGamma
protected

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