Geant4  10.02.p03
G4DNAMillerGreenExcitationModel Class Reference

#include <G4DNAMillerGreenExcitationModel.hh>

Inheritance diagram for G4DNAMillerGreenExcitationModel:
Collaboration 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 G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 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 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 *> *)
 
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=0)
 
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 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

G4ParticleChangeForGamma * fParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4double PartialCrossSection (G4double energy, G4int level, const G4ParticleDefinition *particle)
 
G4double Sum (G4double energy, const G4ParticleDefinition *particle)
 
G4int RandomSelect (G4double energy, const G4ParticleDefinition *particle)
 
G4double S_1s (G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
 
G4double S_2s (G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
 
G4double S_2p (G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
 
G4double R (G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
 
G4DNAMillerGreenExcitationModeloperator= (const G4DNAMillerGreenExcitationModel &right)
 
 G4DNAMillerGreenExcitationModel (const G4DNAMillerGreenExcitationModel &)
 

Private Attributes

G4bool statCode
 
const std::vector< G4double > * fpMolWaterDensity
 
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
 
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
 
G4bool isInitialised
 
G4int verboseLevel
 
G4int nLevels
 
G4DNAWaterExcitationStructure waterExcitation
 
G4double kineticEnergyCorrection [4]
 
G4double slaterEffectiveCharge [3][4]
 
G4double sCoefficient [3][4]
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
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

Definition at line 44 of file G4DNAMillerGreenExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMillerGreenExcitationModel() [1/2]

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");
47 
48  nLevels=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 
68  // Selection of stationary mode
69 
70  statCode = false;
71 }
const std::vector< G4double > * fpMolWaterDensity
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61

◆ ~G4DNAMillerGreenExcitationModel()

G4DNAMillerGreenExcitationModel::~G4DNAMillerGreenExcitationModel ( )
virtual

Definition at line 75 of file G4DNAMillerGreenExcitationModel.cc.

76 {}

◆ G4DNAMillerGreenExcitationModel() [2/2]

G4DNAMillerGreenExcitationModel::G4DNAMillerGreenExcitationModel ( const G4DNAMillerGreenExcitationModel )
private

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 226 of file G4DNAMillerGreenExcitationModel.cc.

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

◆ GetPartialCrossSection()

G4double G4DNAMillerGreenExcitationModel::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particleDefinition,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 407 of file G4DNAMillerGreenExcitationModel.cc.

411 {
412  return PartialCrossSection(kineticEnergy, level, particleDefinition);
413 }
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
Here is the call graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 80 of file G4DNAMillerGreenExcitationModel.cc.

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

◆ operator=()

G4DNAMillerGreenExcitationModel& G4DNAMillerGreenExcitationModel::operator= ( const G4DNAMillerGreenExcitationModel right)
private

◆ PartialCrossSection()

G4double G4DNAMillerGreenExcitationModel::PartialCrossSection ( G4double  energy,
G4int  level,
const G4ParticleDefinition particle 
)
private

Definition at line 417 of file G4DNAMillerGreenExcitationModel.cc.

419 {
420  // ( ( z * aj ) ^ omegaj ) * ( t - ej ) ^ nu
421  // sigma(t) = zEff^2 * sigma0 * --------------------------------------------
422  // jj ^ ( omegaj + nu ) + t ^ ( omegaj + nu )
423  //
424  // where t is the kinetic energy corrected by Helium mass over proton mass for Helium ions
425  //
426  // zEff is:
427  // 1 for protons
428  // 2 for alpha++
429  // and 2 - c1 S_1s - c2 S_2s - c3 S_2p for alpha+ and He
430  //
431  // Dingfelder et al., RPC 59, 255-275, 2000 from Miller and Green (1973)
432  // Formula (34) and Table 2
433 
434  const G4double sigma0(1.E+8 * barn);
435  const G4double nu(1.);
436  const G4double aj[]={876.*eV, 2084.* eV, 1373.*eV, 692.*eV, 900.*eV};
437  const G4double jj[]={19820.*eV, 23490.*eV, 27770.*eV, 30830.*eV, 33080.*eV};
438  const G4double omegaj[]={0.85, 0.88, 0.88, 0.78, 0.78};
439 
440  // Dingfelder's excitation levels
441  const G4double Eliq[5]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
442 
443  G4int particleTypeIndex = 0;
446 
447  if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0;
448  if (particleDefinition == instance->GetIon("hydrogen")) particleTypeIndex=0;
449  if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1;
450  if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2;
451  if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=3;
452 
453  G4double tCorrected;
454  tCorrected = k * kineticEnergyCorrection[particleTypeIndex];
455 
456  // SI - added protection
457  if (tCorrected < Eliq[excitationLevel]) return 0;
458  //
459 
460  G4int z = 10;
461 
462  G4double numerator;
463  numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) *
464  std::pow(tCorrected - Eliq[excitationLevel], nu);
465 
466  // H case : see S. Uehara et al. IJRB 77, 2, 139-154 (2001) - section 3.3
467 
468  if (particleDefinition == instance->GetIon("hydrogen"))
469  numerator = std::pow(z * 0.75*aj[excitationLevel], omegaj[excitationLevel]) *
470  std::pow(tCorrected - Eliq[excitationLevel], nu);
471 
472 
473  G4double power;
474  power = omegaj[excitationLevel] + nu;
475 
476  G4double denominator;
477  denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
478 
479  G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
480 
481  zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, Eliq[excitationLevel], slaterEffectiveCharge[0][particleTypeIndex], 1.) +
482  sCoefficient[1][particleTypeIndex] * S_2s(k, Eliq[excitationLevel], slaterEffectiveCharge[1][particleTypeIndex], 2.) +
483  sCoefficient[2][particleTypeIndex] * S_2p(k, Eliq[excitationLevel], slaterEffectiveCharge[2][particleTypeIndex], 2.) );
484 
485  if (particleDefinition == instance->GetIon("hydrogen")) zEff = 1.;
486 
487  G4double cross = sigma0 * zEff * zEff * numerator / denominator;
488 
489 
490  return cross;
491 }
G4double S_1s(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
static G4DNAGenericIonsManager * Instance(void)
static const double eV
Definition: G4SIunits.hh:212
G4double S_2s(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
static MCTruthManager * instance
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
G4double S_2p(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
G4ParticleDefinition * GetIon(const G4String &name)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ R()

G4double G4DNAMillerGreenExcitationModel::R ( G4double  t,
G4double  energyTransferred,
G4double  slaterEffectiveCharge,
G4double  shellNumber 
)
private

Definition at line 642 of file G4DNAMillerGreenExcitationModel.cc.

646 {
647  // tElectron = m_electron / m_alpha * t
648  // Dingfelder, in Chattanooga 2005 proceedings, p 4
649 
650  G4double tElectron = 0.511/3728. * t;
651 
652  // The following is provided by M. Dingfelder
653  G4double H = 2.*13.60569172 * eV;
654  G4double value = std::sqrt ( 2. * tElectron / H ) / ( energyTransferred / H ) * (_slaterEffectiveCharge/shellNumber);
655 
656  return value;
657 }
static const double eV
Definition: G4SIunits.hh:212
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ RandomSelect()

G4int G4DNAMillerGreenExcitationModel::RandomSelect ( G4double  energy,
const G4ParticleDefinition particle 
)
private

Definition at line 495 of file G4DNAMillerGreenExcitationModel.cc.

496 {
497  G4int i = nLevels;
498  G4double value = 0.;
499  std::deque<double> values;
500 
503 
504  if ( particle == instance->GetIon("alpha++") ||
505  particle == G4Proton::ProtonDefinition()||
506  particle == instance->GetIon("hydrogen") ||
507  particle == instance->GetIon("alpha+") ||
508  particle == instance->GetIon("helium")
509  )
510  {
511  while (i > 0)
512  {
513  i--;
514  G4double partial = PartialCrossSection(k,i,particle);
515  values.push_front(partial);
516  value += partial;
517  }
518 
519  value *= G4UniformRand();
520 
521  i = nLevels;
522 
523  while (i > 0)
524  {
525  i--;
526  if (values[i] > value) return i;
527  value -= values[i];
528  }
529  }
530 
531  /*
532  // add ONE or TWO electron-water excitation for alpha+ and helium
533 
534  if ( particle == instance->GetIon("alpha+")
535  ||
536  particle == instance->GetIon("helium")
537  )
538  {
539  while (i>0)
540  {
541  i--;
542 
543  G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
544  excitationXS->Initialise(G4Electron::ElectronDefinition());
545 
546  G4double sigmaExcitation=0;
547 
548  if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation = excitationXS->PartialCrossSection(k*0.511/3728,i);
549 
550  G4double partial = PartialCrossSection(k,i,particle);
551 
552  if (particle == instance->GetIon("alpha+")) partial = PartialCrossSection(k,i,particle) + sigmaExcitation;
553  if (particle == instance->GetIon("helium")) partial = PartialCrossSection(k,i,particle) + 2*sigmaExcitation;
554 
555  values.push_front(partial);
556  value += partial;
557  delete excitationXS;
558  }
559 
560  value*=G4UniformRand();
561 
562  i=5;
563  while (i>0)
564  {
565  i--;
566 
567  if (values[i]>value) return i;
568 
569  value-=values[i];
570  }
571  }
572 */
573 
574  return 0;
575 }
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
static G4DNAGenericIonsManager * Instance(void)
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
static MCTruthManager * instance
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetIon(const G4String &name)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ S_1s()

G4double G4DNAMillerGreenExcitationModel::S_1s ( G4double  t,
G4double  energyTransferred,
G4double  slaterEffectiveCharge,
G4double  shellNumber 
)
private

Definition at line 592 of file G4DNAMillerGreenExcitationModel.cc.

596 {
597  // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
598  // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
599 
600  G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
601  G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
602 
603  return value;
604 }
G4double R(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ S_2p()

G4double G4DNAMillerGreenExcitationModel::S_2p ( G4double  t,
G4double  energyTransferred,
G4double  slaterEffectiveCharge,
G4double  shellNumber 
)
private

Definition at line 626 of file G4DNAMillerGreenExcitationModel.cc.

630 {
631  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
632  // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
633 
634  G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
635  G4double value = 1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
636 
637  return value;
638 }
G4double R(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ S_2s()

G4double G4DNAMillerGreenExcitationModel::S_2s ( G4double  t,
G4double  energyTransferred,
G4double  slaterEffectiveCharge,
G4double  shellNumber 
)
private

Definition at line 609 of file G4DNAMillerGreenExcitationModel.cc.

613 {
614  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
615  // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
616 
617  G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
618  G4double value = 1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
619 
620  return value;
621 
622 }
G4double R(G4double t, G4double energyTransferred, G4double slaterEffectiveCharge, G4double shellNumber)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 364 of file G4DNAMillerGreenExcitationModel.cc.

369 {
370 
371  if (verboseLevel > 3)
372  G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
373 
374  G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
375 
376  G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
377 
378  // G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
379 
380  // Dingfelder's excitation levels
381  const G4double excitation[]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
382  G4double excitationEnergy = excitation[level];
383 
384  G4double newEnergy = 0.;
385 
386  if (!statCode) newEnergy = particleEnergy0 - excitationEnergy;
387 
388  else newEnergy = particleEnergy0;
389 
390  if (newEnergy>0)
391  {
392  fParticleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
393  fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
394  fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
395 
396  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
398  level,
399  theIncomingTrack);
400 
401  }
402 
403 }
G4int RandomSelect(G4double energy, const G4ParticleDefinition *particle)
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
const G4ThreeVector & GetMomentumDirection() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SelectStationary()

void G4DNAMillerGreenExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 139 of file G4DNAMillerGreenExcitationModel.hh.

140 {
141  statCode = input;
142 }

◆ Sum()

G4double G4DNAMillerGreenExcitationModel::Sum ( G4double  energy,
const G4ParticleDefinition particle 
)
private

Definition at line 579 of file G4DNAMillerGreenExcitationModel.cc.

580 {
581  G4double totalCrossSection = 0.;
582 
583  for (G4int i=0; i<nLevels; i++)
584  {
585  totalCrossSection += PartialCrossSection(k,i,particle);
586  }
587  return totalCrossSection;
588 }
int G4int
Definition: G4Types.hh:78
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMillerGreenExcitationModel::fParticleChangeForGamma
protected

Definition at line 77 of file G4DNAMillerGreenExcitationModel.hh.

◆ fpMolWaterDensity

const std::vector<G4double>* G4DNAMillerGreenExcitationModel::fpMolWaterDensity
private

Definition at line 84 of file G4DNAMillerGreenExcitationModel.hh.

◆ highEnergyLimit

std::map<G4String,G4double,std::less<G4String> > G4DNAMillerGreenExcitationModel::highEnergyLimit
private

Definition at line 87 of file G4DNAMillerGreenExcitationModel.hh.

◆ isInitialised

G4bool G4DNAMillerGreenExcitationModel::isInitialised
private

Definition at line 89 of file G4DNAMillerGreenExcitationModel.hh.

◆ kineticEnergyCorrection

G4double G4DNAMillerGreenExcitationModel::kineticEnergyCorrection[4]
private

Definition at line 126 of file G4DNAMillerGreenExcitationModel.hh.

◆ lowEnergyLimit

std::map<G4String,G4double,std::less<G4String> > G4DNAMillerGreenExcitationModel::lowEnergyLimit
private

Definition at line 86 of file G4DNAMillerGreenExcitationModel.hh.

◆ nLevels

G4int G4DNAMillerGreenExcitationModel::nLevels
private

Definition at line 102 of file G4DNAMillerGreenExcitationModel.hh.

◆ sCoefficient

G4double G4DNAMillerGreenExcitationModel::sCoefficient[3][4]
private

Definition at line 128 of file G4DNAMillerGreenExcitationModel.hh.

◆ slaterEffectiveCharge

G4double G4DNAMillerGreenExcitationModel::slaterEffectiveCharge[3][4]
private

Definition at line 127 of file G4DNAMillerGreenExcitationModel.hh.

◆ statCode

G4bool G4DNAMillerGreenExcitationModel::statCode
private

Definition at line 81 of file G4DNAMillerGreenExcitationModel.hh.

◆ verboseLevel

G4int G4DNAMillerGreenExcitationModel::verboseLevel
private

Definition at line 90 of file G4DNAMillerGreenExcitationModel.hh.

◆ waterExcitation

G4DNAWaterExcitationStructure G4DNAMillerGreenExcitationModel::waterExcitation
private

Definition at line 104 of file G4DNAMillerGreenExcitationModel.hh.


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