Geant4  10.02.p03
G4DNADingfelderChargeDecreaseModel Class Reference

#include <G4DNADingfelderChargeDecreaseModel.hh>

Inheritance diagram for G4DNADingfelderChargeDecreaseModel:
Collaboration diagram for G4DNADingfelderChargeDecreaseModel:

Public Member Functions

 G4DNADingfelderChargeDecreaseModel (const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeDecreaseModel")
 
virtual ~G4DNADingfelderChargeDecreaseModel ()
 
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 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 GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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)
 
G4int NumberOfFinalStates (G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
 
G4ParticleDefinitionOutgoingParticleDefinition (G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
 
G4double WaterBindingEnergyConstant (G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
 
G4double OutgoingParticleBindingEnergyConstant (G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
 
G4DNADingfelderChargeDecreaseModeloperator= (const G4DNADingfelderChargeDecreaseModel &right)
 
 G4DNADingfelderChargeDecreaseModel (const G4DNADingfelderChargeDecreaseModel &)
 

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 numberOfPartialCrossSections [3]
 
G4double f0 [2][3]
 
G4double a0 [2][3]
 
G4double a1 [2][3]
 
G4double b0 [2][3]
 
G4double b1 [2][3]
 
G4double c0 [2][3]
 
G4double d0 [2][3]
 
G4double x0 [2][3]
 
G4double x1 [2][3]
 

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 40 of file G4DNADingfelderChargeDecreaseModel.hh.

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeDecreaseModel() [1/2]

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeDecreaseModel" 
)

Definition at line 41 of file G4DNADingfelderChargeDecreaseModel.cc.

42  :
43  G4VEmModel(nam), isInitialised(false)
44 {
45  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
50 
51  verboseLevel = 0;
52  // Verbosity scale:
53  // 0 = nothing
54  // 1 = warning for energy non-conservation
55  // 2 = details of energy budget
56  // 3 = calculation of cross sections, file openings, sampling of atoms
57  // 4 = entering in methods
58 
59  if (verboseLevel > 0)
60  {
61  G4cout << "Dingfelder charge decrease model is constructed " << G4endl;
62  }
64 
65  // Selection of stationary mode
66 
67  statCode = false;
68 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

◆ ~G4DNADingfelderChargeDecreaseModel()

G4DNADingfelderChargeDecreaseModel::~G4DNADingfelderChargeDecreaseModel ( )
virtual

Definition at line 72 of file G4DNADingfelderChargeDecreaseModel.cc.

73 {
74 }

◆ G4DNADingfelderChargeDecreaseModel() [2/2]

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4DNADingfelderChargeDecreaseModel )
private

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 211 of file G4DNADingfelderChargeDecreaseModel.cc.

216 {
217  if (verboseLevel > 3)
218  {
219  G4cout
220  << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeDecreaseModel"
221  << G4endl;
222  }
223 
224  // Calculate total cross section for model
225 
228 
229  if (
230  particleDefinition != G4Proton::ProtonDefinition()
231  &&
232  particleDefinition != instance->GetIon("alpha++")
233  &&
234  particleDefinition != instance->GetIon("alpha+")
235  )
236 
237  return 0;
238 
239  G4double lowLim = 0;
240  G4double highLim = 0;
241  G4double crossSection = 0.;
242 
243  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
244 
245  if(waterDensity!= 0.0)
246  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
247  {
248  const G4String& particleName = particleDefinition->GetParticleName();
249 
250  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
251  pos1 = lowEnergyLimit.find(particleName);
252 
253  if (pos1 != lowEnergyLimit.end())
254  {
255  lowLim = pos1->second;
256  }
257 
258  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
259  pos2 = highEnergyLimit.find(particleName);
260 
261  if (pos2 != highEnergyLimit.end())
262  {
263  highLim = pos2->second;
264  }
265 
266  if (k >= lowLim && k < highLim)
267  {
268  crossSection = Sum(k,particleDefinition);
269  }
270 
271  if (verboseLevel > 2)
272  {
273  G4cout << "_______________________________________" << G4endl;
274  G4cout << "G4DNADingfelderChargeDecreaeModel" << G4endl;
275  G4cout << "Kinetic energy(eV)=" << k/eV << "particle :" << particleName << G4endl;
276  G4cout << "Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
277  G4cout << "Cross section per water molecule (cm^-1)=" << crossSection*
278  waterDensity/(1./cm) << G4endl;
279 // material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
280  }
281  }
282 
283  return crossSection*waterDensity;
284 // return crossSection*material->GetAtomicNumDensityVector()[1];
285 
286 }
static const double cm
Definition: G4SIunits.hh:118
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4GLOB_DLL std::ostream G4cout
static G4DNAGenericIonsManager * Instance(void)
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
static const double eV
Definition: G4SIunits.hh:212
G4double Sum(G4double energy, const G4ParticleDefinition *particle)
#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:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 78 of file G4DNADingfelderChargeDecreaseModel.cc.

80 {
81 
82  if (verboseLevel > 3)
83  {
84  G4cout << "Calling G4DNADingfelderChargeDecreaseModel::Initialise()"
85  << G4endl;
86  }
87 
88  // Energy limits
89 
93  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
94  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
95 
97  G4String alphaPlusPlus;
98  G4String alphaPlus;
99 
100  // LIMITS
101 
102  proton = protonDef->GetParticleName();
103  lowEnergyLimit[proton] = 100. * eV;
104  highEnergyLimit[proton] = 100. * MeV;
105 
106  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
107  lowEnergyLimit[alphaPlusPlus] = 1. * keV;
108  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
109 
110  alphaPlus = alphaPlusDef->GetParticleName();
111  lowEnergyLimit[alphaPlus] = 1. * keV;
112  highEnergyLimit[alphaPlus] = 400. * MeV;
113 
114  //
115 
116  if (particle==protonDef)
117  {
120  }
121 
122  if (particle==alphaPlusPlusDef)
123  {
124  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
125  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
126  }
127 
128  if (particle==alphaPlusDef)
129  {
130  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
132  }
133 
134  // Final state
135 
136  //PROTON
137  f0[0][0]=1.;
138  a0[0][0]=-0.180;
139  a1[0][0]=-3.600;
140  b0[0][0]=-18.22;
141  b1[0][0]=-1.997;
142  c0[0][0]=0.215;
143  d0[0][0]=3.550;
144  x0[0][0]=3.450;
145  x1[0][0]=5.251;
146 
148 
149  //ALPHA++
150  f0[0][1]=1.; a0[0][1]=0.95;
151  a1[0][1]=-2.75;
152  b0[0][1]=-23.00;
153  c0[0][1]=0.215;
154  d0[0][1]=2.95;
155  x0[0][1]=3.50;
156 
157  f0[1][1]=1.;
158  a0[1][1]=0.95;
159  a1[1][1]=-2.75;
160  b0[1][1]=-23.73;
161  c0[1][1]=0.250;
162  d0[1][1]=3.55;
163  x0[1][1]=3.72;
164 
165  x1[0][1]=-1.;
166  b1[0][1]=-1.;
167 
168  x1[1][1]=-1.;
169  b1[1][1]=-1.;
170 
172 
173  // ALPHA+
174  f0[0][2]=1.;
175  a0[0][2]=0.65;
176  a1[0][2]=-2.75;
177  b0[0][2]=-21.81;
178  c0[0][2]=0.232;
179  d0[0][2]=2.95;
180  x0[0][2]=3.53;
181 
182  x1[0][2]=-1.;
183  b1[0][2]=-1.;
184 
186 
187  //
188 
189  if( verboseLevel>0 )
190  {
191  G4cout << "Dingfelder charge decrease model is initialized " << G4endl
192  << "Energy range: "
193  << LowEnergyLimit() / keV << " keV - "
194  << HighEnergyLimit() / MeV << " MeV for "
195  << particle->GetParticleName()
196  << G4endl;
197  }
198 
199  // Initialize water density pointer
201 
202  if (isInitialised)
203  { return;}
205  isInitialised = true;
206 
207 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static G4DNAGenericIonsManager * Instance(void)
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
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:

◆ NumberOfFinalStates()

G4int G4DNADingfelderChargeDecreaseModel::NumberOfFinalStates ( G4ParticleDefinition particleDefinition,
G4int  finalStateIndex 
)
private

Definition at line 367 of file G4DNADingfelderChargeDecreaseModel.cc.

370 {
371  if (particleDefinition == G4Proton::Proton())
372  return 1;
373 
376 
377  if (particleDefinition == instance->GetIon("alpha++"))
378  {
379  if (finalStateIndex == 0)
380  return 1;
381  return 2;
382  }
383 
384  if (particleDefinition == instance->GetIon("alpha+"))
385  return 1;
386 
387  return 0;
388 }
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
static MCTruthManager * instance
G4ParticleDefinition * GetIon(const G4String &name)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ OutgoingParticleBindingEnergyConstant()

G4double G4DNADingfelderChargeDecreaseModel::OutgoingParticleBindingEnergyConstant ( G4ParticleDefinition particleDefinition,
G4int  finalStateIndex 
)
private

Definition at line 457 of file G4DNADingfelderChargeDecreaseModel.cc.

459 {
461 
462  if (particleDefinition == G4Proton::Proton())
463  return 13.6 * eV;
464 
465  if (particleDefinition == instance->GetIon("alpha++"))
466  {
467  // Binding energy for He+ -> He++ + e- 54.509 eV
468  // Binding energy for He -> He+ + e- 24.587 eV
469 
470  if (finalStateIndex == 0)
471  return 54.509 * eV;
472 
473  return (54.509 + 24.587) * eV;
474  }
475 
476  if (particleDefinition == instance->GetIon("alpha+"))
477  {
478  // Binding energy for He+ -> He++ + e- 54.509 eV
479  // Binding energy for He -> He+ + e- 24.587 eV
480 
481  return 24.587 * eV;
482  }
483 
484  return 0;
485 }
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
static const double eV
Definition: G4SIunits.hh:212
static MCTruthManager * instance
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OutgoingParticleDefinition()

G4ParticleDefinition * G4DNADingfelderChargeDecreaseModel::OutgoingParticleDefinition ( G4ParticleDefinition particleDefinition,
G4int  finalStateIndex 
)
private

Definition at line 392 of file G4DNADingfelderChargeDecreaseModel.cc.

394 {
396 
397  if (particleDefinition == G4Proton::Proton())
398  return instance->GetIon("hydrogen");
399 
400  if (particleDefinition == instance->GetIon("alpha++"))
401  {
402  if (finalStateIndex == 0)
403  return instance->GetIon("alpha+");
404  return instance->GetIon("helium");
405  }
406 
407  if (particleDefinition == instance->GetIon("alpha+"))
408  return instance->GetIon("helium");
409 
410  return 0;
411 }
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
static MCTruthManager * instance
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PartialCrossSection()

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

Definition at line 489 of file G4DNADingfelderChargeDecreaseModel.cc.

492 {
493  G4int particleTypeIndex = 0;
496 
497  if (particleDefinition == G4Proton::ProtonDefinition())
498  particleTypeIndex = 0;
499  if (particleDefinition == instance->GetIon("alpha++"))
500  particleTypeIndex = 1;
501  if (particleDefinition == instance->GetIon("alpha+"))
502  particleTypeIndex = 2;
503 
504  //
505  // sigma(T) = f0 10 ^ y(log10(T/eV))
506  //
507  // / a0 x + b0 x < x0
508  // |
509  // y(x) = < a0 x + b0 - c0 (x - x0)^d0 x0 <= x < x1
510  // |
511  // \ a1 x + b1 x >= x1
512  //
513  //
514  // f0, a0, a1, b0, b1, c0, d0, x0, x1 are parameters that change for protons and helium (0, +, ++)
515  //
516  // f0 has been added to the code in order to manage partial (shell-dependent) cross sections (if no shell dependence is present. f0=1. Sum of f0 over the considered shells should give 1)
517  //
518  // From Rad. Phys. and Chem. 59 (2000) 255-275, M. Dingfelder et al.
519  // Inelastic-collision cross sections of liquid water for interactions of energetic proton
520  //
521 
522  if (x1[index][particleTypeIndex] < x0[index][particleTypeIndex])
523  {
524  //
525  // if x1 < x0 means that x1 and b1 will be calculated with the following formula (this piece of code is run on all alphas and not on protons)
526  //
527  // x1 = x0 + ((a0 - a1)/(c0 * d0)) ^ (1 / (d0 - 1))
528  //
529  // b1 = (a0 - a1) * x1 + b0 - c0 * (x1 - x0) ^ d0
530  //
531 
532  x1[index][particleTypeIndex] = x0[index][particleTypeIndex]
533  + std::pow((a0[index][particleTypeIndex] - a1[index][particleTypeIndex])
534  / (c0[index][particleTypeIndex]
535  * d0[index][particleTypeIndex]),
536  1. / (d0[index][particleTypeIndex] - 1.));
537  b1[index][particleTypeIndex] = (a0[index][particleTypeIndex]
538  - a1[index][particleTypeIndex]) * x1[index][particleTypeIndex]
539  + b0[index][particleTypeIndex]
540  - c0[index][particleTypeIndex]
541  * std::pow(x1[index][particleTypeIndex]
542  - x0[index][particleTypeIndex],
543  d0[index][particleTypeIndex]);
544  }
545 
546  G4double x(std::log10(k / eV));
547  G4double y;
548 
549  if (x < x0[index][particleTypeIndex])
550  y = a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex];
551  else if (x < x1[index][particleTypeIndex])
552  y = a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex]
553  - c0[index][particleTypeIndex]
554  * std::pow(x - x0[index][particleTypeIndex],
555  d0[index][particleTypeIndex]);
556  else
557  y = a1[index][particleTypeIndex] * x + b1[index][particleTypeIndex];
558 
559  return f0[index][particleTypeIndex] * std::pow(10., y) * m * m;
560 
561 }
Int_t index
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
Double_t y
static G4DNAGenericIonsManager * Instance(void)
static const double eV
Definition: G4SIunits.hh:212
static const double m
Definition: G4SIunits.hh:128
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:

◆ RandomSelect()

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

Definition at line 563 of file G4DNADingfelderChargeDecreaseModel.cc.

565 {
566  G4int particleTypeIndex = 0;
569 
570  if (particleDefinition == G4Proton::ProtonDefinition())
571  particleTypeIndex = 0;
572  if (particleDefinition == instance->GetIon("alpha++"))
573  particleTypeIndex = 1;
574  if (particleDefinition == instance->GetIon("alpha+"))
575  particleTypeIndex = 2;
576 
577  const G4int n = numberOfPartialCrossSections[particleTypeIndex];
578  G4double* values(new G4double[n]);
579  G4double value(0);
580  G4int i = n;
581 
582  while (i > 0)
583  {
584  i--;
585  values[i] = PartialCrossSection(k, i, particleDefinition);
586  value += values[i];
587  }
588 
589  value *= G4UniformRand();
590 
591  i = n;
592  while (i > 0)
593  {
594  i--;
595 
596  if (values[i] > value)
597  break;
598 
599  value -= values[i];
600  }
601 
602  delete[] values;
603 
604  return i;
605 }
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
Char_t n[5]
#define G4UniformRand()
Definition: Randomize.hh:97
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
static G4DNAGenericIonsManager * Instance(void)
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:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 290 of file G4DNADingfelderChargeDecreaseModel.cc.

296 {
297  if (verboseLevel > 3)
298  {
299  G4cout
300  << "Calling SampleSecondaries() of G4DNADingfelderChargeDecreaseModel"
301  << G4endl;
302  }
303 
304  G4double inK = aDynamicParticle->GetKineticEnergy();
305 
306  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
307 
308  G4double particleMass = definition->GetPDGMass();
309 
310  G4int finalStateIndex = RandomSelect(inK,definition);
311 
312  G4int n = NumberOfFinalStates(definition, finalStateIndex);
313  G4double waterBindingEnergy = WaterBindingEnergyConstant(definition, finalStateIndex);
314  G4double outgoingParticleBindingEnergy = OutgoingParticleBindingEnergyConstant(definition, finalStateIndex);
315 
316  G4double outK = 0.;
317 
318  if (definition==G4Proton::Proton())
319  {
320  if (!statCode) outK = inK - n*(inK*electron_mass_c2/proton_mass_c2)
321  - waterBindingEnergy + outgoingParticleBindingEnergy;
322  else outK = inK;
323  }
324 
325  else
326  {
327  if (!statCode) outK = inK - n*(inK*electron_mass_c2/particleMass)
328  - waterBindingEnergy + outgoingParticleBindingEnergy;
329  else outK = inK;
330  }
331 
332  if (outK<0)
333  {
334  G4Exception("G4DNADingfelderChargeDecreaseModel::SampleSecondaries","em0004",
335  FatalException,"Final kinetic energy is negative.");
336  }
337 
338  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
339 
340  if (!statCode) fParticleChangeForGamma->ProposeLocalEnergyDeposit(waterBindingEnergy);
341 
342  else
343 
344  {
345  if (definition==G4Proton::Proton())
346  fParticleChangeForGamma->ProposeLocalEnergyDeposit(n*(inK*electron_mass_c2/proton_mass_c2)
347  + waterBindingEnergy - outgoingParticleBindingEnergy);
348  else
349  fParticleChangeForGamma->ProposeLocalEnergyDeposit(n*(inK*electron_mass_c2/particleMass)
350  + waterBindingEnergy - outgoingParticleBindingEnergy);
351  }
352 
353  G4DynamicParticle* dp = new G4DynamicParticle (OutgoingParticleDefinition(definition, finalStateIndex),
354  aDynamicParticle->GetMomentumDirection(),
355  outK);
356  fvect->push_back(dp);
357 
358  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
360  1,
361  theIncomingTrack);
362 
363 }
int G4int
Definition: G4Types.hh:78
G4double WaterBindingEnergyConstant(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
Char_t n[5]
G4double GetKineticEnergy() const
G4double OutgoingParticleBindingEnergyConstant(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
const G4ThreeVector & GetMomentumDirection() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4int NumberOfFinalStates(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4ParticleDefinition * OutgoingParticleDefinition(G4ParticleDefinition *particleDefinition, G4int finalStateIndex)
double G4double
Definition: G4Types.hh:76
G4int RandomSelect(G4double energy, const G4ParticleDefinition *particle)
Here is the call graph for this function:

◆ SelectStationary()

void G4DNADingfelderChargeDecreaseModel::SelectStationary ( G4bool  input)
inline

Definition at line 122 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ Sum()

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

Definition at line 609 of file G4DNADingfelderChargeDecreaseModel.cc.

611 {
612  G4int particleTypeIndex = 0;
615 
616  if (particleDefinition == G4Proton::ProtonDefinition())
617  particleTypeIndex = 0;
618  if (particleDefinition == instance->GetIon("alpha++"))
619  particleTypeIndex = 1;
620  if (particleDefinition == instance->GetIon("alpha+"))
621  particleTypeIndex = 2;
622 
623  G4double totalCrossSection = 0.;
624 
625  for (G4int i = 0; i < numberOfPartialCrossSections[particleTypeIndex]; i++)
626  {
627  totalCrossSection += PartialCrossSection(k, i, particleDefinition);
628  }
629  return totalCrossSection;
630 }
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition *particle)
static G4DNAGenericIonsManager * Instance(void)
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:

◆ WaterBindingEnergyConstant()

G4double G4DNADingfelderChargeDecreaseModel::WaterBindingEnergyConstant ( G4ParticleDefinition particleDefinition,
G4int  finalStateIndex 
)
private

Definition at line 415 of file G4DNADingfelderChargeDecreaseModel.cc.

417 {
418  // Ionization energy of first water shell
419  // Rad. Phys. Chem. 59 p.267 by Dingf. et al.
420  // W + 10.79 eV -> W+ + e-
421 
423 
424  if (particleDefinition == G4Proton::Proton())
425  return 10.79 * eV;
426 
427  if (particleDefinition == instance->GetIon("alpha++"))
428  {
429  // Binding energy for W+ -> W++ + e- 10.79 eV
430  // Binding energy for W -> W+ + e- 10.79 eV
431  //
432  // Ionization energy of first water shell
433  // Rad. Phys. Chem. 59 p.267 by Dingf. et al.
434 
435  if (finalStateIndex == 0)
436  return 10.79 * eV;
437 
438  return 10.79 * 2 * eV;
439  }
440 
441  if (particleDefinition == instance->GetIon("alpha+"))
442  {
443  // Binding energy for W+ -> W++ + e- 10.79 eV
444  // Binding energy for W -> W+ + e- 10.79 eV
445  //
446  // Ionization energy of first water shell
447  // Rad. Phys. Chem. 59 p.267 by Dingf. et al.
448 
449  return 10.79 * eV;
450  }
451 
452  return 0;
453 }
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
static const double eV
Definition: G4SIunits.hh:212
static MCTruthManager * instance
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ a0

G4double G4DNADingfelderChargeDecreaseModel::a0[2][3]
private

Definition at line 94 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ a1

G4double G4DNADingfelderChargeDecreaseModel::a1[2][3]
private

Definition at line 95 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ b0

G4double G4DNADingfelderChargeDecreaseModel::b0[2][3]
private

Definition at line 96 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ b1

G4double G4DNADingfelderChargeDecreaseModel::b1[2][3]
private

Definition at line 97 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ c0

G4double G4DNADingfelderChargeDecreaseModel::c0[2][3]
private

Definition at line 98 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ d0

G4double G4DNADingfelderChargeDecreaseModel::d0[2][3]
private

Definition at line 99 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ f0

G4double G4DNADingfelderChargeDecreaseModel::f0[2][3]
private

Definition at line 93 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeDecreaseModel::fParticleChangeForGamma
protected

Definition at line 68 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ fpMolWaterDensity

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

Definition at line 75 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ highEnergyLimit

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

Definition at line 78 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ isInitialised

G4bool G4DNADingfelderChargeDecreaseModel::isInitialised
private

Definition at line 80 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ lowEnergyLimit

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

Definition at line 77 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ numberOfPartialCrossSections

G4int G4DNADingfelderChargeDecreaseModel::numberOfPartialCrossSections[3]
private

Definition at line 91 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ statCode

G4bool G4DNADingfelderChargeDecreaseModel::statCode
private

Definition at line 72 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ verboseLevel

G4int G4DNADingfelderChargeDecreaseModel::verboseLevel
private

Definition at line 81 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ x0

G4double G4DNADingfelderChargeDecreaseModel::x0[2][3]
private

Definition at line 100 of file G4DNADingfelderChargeDecreaseModel.hh.

◆ x1

G4double G4DNADingfelderChargeDecreaseModel::x1[2][3]
private

Definition at line 101 of file G4DNADingfelderChargeDecreaseModel.hh.


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