Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAEmfietzoglouIonisationModel Class Reference

#include <G4DNAEmfietzoglouIonisationModel.hh>

Inheritance diagram for G4DNAEmfietzoglouIonisationModel:
Collaboration diagram for G4DNAEmfietzoglouIonisationModel:

Public Member Functions

 G4DNAEmfietzoglouIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
 
virtual ~G4DNAEmfietzoglouIonisationModel ()
 
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)
 
double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
void SelectFasterComputation (G4bool input)
 
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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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=nullptr)
 
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

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 54 of file G4DNAEmfietzoglouIonisationModel.hh.

Constructor & Destructor Documentation

G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAEmfietzoglouIonisationModel" 
)

Definition at line 52 of file G4DNAEmfietzoglouIonisationModel.cc.

53  :
54  G4VEmModel(nam), isInitialised(false)
55 {
56  verboseLevel = 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, sampling of atoms
62  // 4 = entering in methods
63 
64  if(verboseLevel > 0)
65  {
66  G4cout << "Emfietzoglou ionisation model is constructed " << G4endl;
67  }
68 
69  //Mark this model as "applicable" for atomic deexcitation
70  SetDeexcitationFlag(true);
71  fAtomDeexcitation = 0;
73  fpMolWaterDensity = 0;
74 
75  // define default angular generator
77 
78  SetLowEnergyLimit(10. * eV);
79  SetHighEnergyLimit(10. * keV);
80 
81  // Selection of computation method
82 
83  fasterCode = false;
84 
85  // Selection of stationary mode
86 
87  statCode = false;
88 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4DNAEmfietzoglouIonisationModel::~G4DNAEmfietzoglouIonisationModel ( )
virtual

Definition at line 92 of file G4DNAEmfietzoglouIonisationModel.cc.

93 {
94  // Cross section
95 
96  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
97  for(pos = tableData.begin(); pos != tableData.end(); ++pos)
98  {
99  G4DNACrossSectionDataSet* table = pos->second;
100  delete table;
101  }
102 
103  // Final state
104 
105  eVecm.clear();
106 
107 }
static const G4double pos

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 246 of file G4DNAEmfietzoglouIonisationModel.cc.

251 {
252  if(verboseLevel > 3)
253  {
254  G4cout
255  << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouIonisationModel"
256  << G4endl;
257  }
258 
259  if (particleDefinition != G4Electron::ElectronDefinition()) return 0; // necessary ??
260 
261  // Calculate total cross section for model
262 
263  G4double lowLim = 0;
264  G4double highLim = 0;
265  G4double sigma=0;
266 
267  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
268 
269  if(waterDensity!= 0.0)
270  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
271  {
272  const G4String& particleName = particleDefinition->GetParticleName();
273 
274  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
275  pos1 = lowEnergyLimit.find(particleName);
276  if (pos1 != lowEnergyLimit.end())
277  {
278  lowLim = pos1->second;
279  }
280 
281  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
282  pos2 = highEnergyLimit.find(particleName);
283  if (pos2 != highEnergyLimit.end())
284  {
285  highLim = pos2->second;
286  }
287 
288  if (ekin >= lowLim && ekin < highLim)
289  {
290  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
291  pos = tableData.find(particleName);
292 
293  if (pos != tableData.end())
294  {
295  G4DNACrossSectionDataSet* table = pos->second;
296  if (table != 0)
297  {
298  sigma = table->FindValue(ekin);
299  }
300  }
301  else
302  {
303  G4Exception("G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume","em0002",
304  FatalException,"Model not applicable to particle type.");
305  }
306  }
307 
308  if (verboseLevel > 2)
309  {
310  G4cout << "__________________________________" << G4endl;
311  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO START" << G4endl;
312  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
313  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
314  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
315  G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO END" << G4endl;
316  }
317  } // if (waterMaterial)
318 
319  return sigma*waterDensity;
320 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos

Here is the call graph for this function:

double G4DNAEmfietzoglouIonisationModel::DifferentialCrossSection ( G4ParticleDefinition aParticleDefinition,
G4double  k,
G4double  energyTransfer,
G4int  shell 
)

Definition at line 594 of file G4DNAEmfietzoglouIonisationModel.cc.

598 {
599  G4double sigma = 0.;
600 
601  if(energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
602  {
603  G4double valueT1 = 0;
604  G4double valueT2 = 0;
605  G4double valueE21 = 0;
606  G4double valueE22 = 0;
607  G4double valueE12 = 0;
608  G4double valueE11 = 0;
609 
610  G4double xs11 = 0;
611  G4double xs12 = 0;
612  G4double xs21 = 0;
613  G4double xs22 = 0;
614 
615  if(particleDefinition == G4Electron::ElectronDefinition())
616  {
617  // k should be in eV and energy transfer eV also
618 
619  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
620  eTdummyVec.end(),
621  k);
622 
623  std::vector<double>::iterator t1 = t2 - 1;
624 
625  // SI : the following condition avoids situations where energyTransfer >last vector element
626  if(energyTransfer <= eVecm[(*t1)].back() && energyTransfer
627  <= eVecm[(*t2)].back())
628  {
629  std::vector<double>::iterator e12 =
630  std::upper_bound(eVecm[(*t1)].begin(),
631  eVecm[(*t1)].end(),
632  energyTransfer);
633  std::vector<double>::iterator e11 = e12 - 1;
634 
635  std::vector<double>::iterator e22 =
636  std::upper_bound(eVecm[(*t2)].begin(),
637  eVecm[(*t2)].end(),
638  energyTransfer);
639  std::vector<double>::iterator e21 = e22 - 1;
640 
641  valueT1 = *t1;
642  valueT2 = *t2;
643  valueE21 = *e21;
644  valueE22 = *e22;
645  valueE12 = *e12;
646  valueE11 = *e11;
647 
648  xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
649  xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
650  xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
651  xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
652 
653  //G4cout << "-------------------" << G4endl;
654  //G4cout << "ionizationLevelIndex=" << ionizationLevelIndex << G4endl;
655  //G4cout << "valueT1/eV=" << valueT1 << " valueT2/eV=" << valueT2 << G4endl;
656  //G4cout << "valueE11/eV=" << valueE11 << " valueE12/eV=" << valueE12 << " valueE21/eV=" << valueE21 << " valueE22/eV=" << valueE22 << G4endl;
657  //G4cout << "xs11=" << xs11 / ((1.e-22 / 3.343) * m*m) << G4endl;
658  //G4cout << "xs12=" << xs12 / ((1.e-22 / 3.343) * m*m) << G4endl;
659  //G4cout << "xs21=" << xs21 / ((1.e-22 / 3.343) * m*m) << G4endl;
660  //G4cout << "xs22=" << xs22 / ((1.e-22 / 3.343) * m*m) << G4endl;
661  //G4cout << "###################" << G4endl;
662 
663  }
664 
665  }
666 
667  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
668  if(xsProduct != 0.)
669  {
670  sigma = QuadInterpolator(valueE11,
671  valueE12,
672  valueE21,
673  valueE22,
674  xs11,
675  xs12,
676  xs21,
677  xs22,
678  valueT1,
679  valueT2,
680  k,
681  energyTransfer);
682  }
683 
684  }
685 
686  return sigma;
687 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 111 of file G4DNAEmfietzoglouIonisationModel.cc.

113 {
114 
115  if(verboseLevel > 3)
116  {
117  G4cout << "Calling G4DNAEmfietzoglouIonisationModel::Initialise()" << G4endl;
118  }
119 
120  // Energy limits
121 
122  G4String fileElectron("dna/sigma_ionisation_e_emfietzoglou");
123 
125 
127 
128  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
129 
130  char *path = getenv("G4LEDATA");
131 
132  // *** ELECTRON
133 
134  electron = electronDef->GetParticleName();
135 
136  tableFile[electron] = fileElectron;
137 
138  lowEnergyLimit[electron] = 10. * eV;
139  highEnergyLimit[electron] = 10. * keV;
140 
141  // Cross section
142 
144  tableE->LoadData(fileElectron);
145 
146  tableData[electron] = tableE;
147 
148  // Final state
149 
150  std::ostringstream eFullFileName;
151 
152  if (fasterCode) eFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_e_emfietzoglou.dat";
153  if (!fasterCode) eFullFileName << path << "/dna/sigmadiff_ionisation_e_emfietzoglou.dat";
154 
155  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
156 
157  if (!eDiffCrossSection)
158  {
159  if (fasterCode) G4Exception("G4DNAEmfietzoglouIonisationModel::Initialise","em0003",
160  FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_e_emfietzoglou.dat");
161 
162  if (!fasterCode) G4Exception("G4DNAEmfietzoglouIonisationModel::Initialise","em0003",
163  FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_emfietzoglou.dat");
164  }
165 
166  //
167 
168  // Clear the arrays for re-initialization case (MT mode)
169  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
170 
171  eTdummyVec.clear();
172 
173  eVecm.clear();
174 
175  eProbaShellMap->clear();
176 
177  eDiffCrossSectionData->clear();
178 
179  eNrjTransfData->clear();
180 
181  //
182 
183  eTdummyVec.push_back(0.);
184  while(!eDiffCrossSection.eof())
185  {
186  double tDummy;
187  double eDummy;
188  eDiffCrossSection>>tDummy>>eDummy;
189  if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
190  for (int j=0; j<5; j++)
191  {
192  eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
193 
194  if (fasterCode)
195  {
196  eNrjTransfData[j][tDummy][eDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
197  eProbaShellMap[j][tDummy].push_back(eDiffCrossSectionData[j][tDummy][eDummy]);
198  }
199 
200  // SI - only if eof is not reached
201  if (!eDiffCrossSection.eof() && !fasterCode)
202  {
203  eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
204  }
205 
206  if (!fasterCode) eVecm[tDummy].push_back(eDummy);
207 
208  }
209  }
210 
211  //
212 
213  if (particle==electronDef)
214  {
215  SetLowEnergyLimit(lowEnergyLimit[electron]);
216  SetHighEnergyLimit(highEnergyLimit[electron]);
217  }
218 
219  if( verboseLevel>0 )
220  {
221  G4cout << "Emfietzoglou ionisation model is initialized " << G4endl
222  << "Energy range: "
223  << LowEnergyLimit() / eV << " eV - "
224  << HighEnergyLimit() / keV << " keV for "
225  << particle->GetParticleName()
226  << G4endl;
227  }
228 
229  // Initialize water density pointer
230  fpMolWaterDensity =
232  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
233 
234  //
235  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
236 
237  if (isInitialised)
238  { return;}
240  isInitialised = true;
241 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
static G4LossTableManager * Instance()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double eV
Definition: G4SIunits.hh:215
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

void G4DNAEmfietzoglouIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 325 of file G4DNAEmfietzoglouIonisationModel.cc.

330 {
331 
332  if(verboseLevel > 3)
333  {
334  G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouIonisationModel"
335  << G4endl;
336  }
337 
338  G4double lowLim = 0;
339  G4double highLim = 0;
340 
341  G4double k = particle->GetKineticEnergy();
342 
343  const G4String& particleName = particle->GetDefinition()->GetParticleName();
344 
345  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
346  pos1 = lowEnergyLimit.find(particleName);
347 
348  if (pos1 != lowEnergyLimit.end())
349  {
350  lowLim = pos1->second;
351  }
352 
353  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
354  pos2 = highEnergyLimit.find(particleName);
355 
356  if (pos2 != highEnergyLimit.end())
357  {
358  highLim = pos2->second;
359  }
360 
361  if (k >= lowLim && k < highLim)
362  {
363  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
364  G4double particleMass = particle->GetDefinition()->GetPDGMass();
365  G4double totalEnergy = k + particleMass;
366  G4double pSquare = k * (totalEnergy + particleMass);
367  G4double totalMomentum = std::sqrt(pSquare);
368 
369  G4int ionizationShell = 0;
370 
371  ionizationShell = RandomSelect(k,particleName);
372 
373  // AM: sample deexcitation
374  // here we assume that H_{2}O electronic levels are the same as Oxygen.
375  // this can be considered true with a rough 10% error in energy on K-shell,
376 
377  G4int secNumberInit = 0;// need to know at a certain point the energy of secondaries
378  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
379 
381  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
382 
383  //SI: additional protection if tcs interpolation method is modified
384  if (k<bindingEnergy) return;
385  //
386 
387  G4int Z = 8;
388  if(fAtomDeexcitation)
389  {
391 
392  if (ionizationShell <5 && ionizationShell >1)
393  {
394  as = G4AtomicShellEnumerator(4-ionizationShell);
395  }
396  else if (ionizationShell <2)
397  {
398  as = G4AtomicShellEnumerator(3);
399  }
400 
401  // FOR DEBUG ONLY
402  // if (ionizationShell == 4) {
403  //
404  // G4cout << "Z: " << Z << " as: " << as
405  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
406  // G4cout << "Press <Enter> key to continue..." << G4endl;
407  // G4cin.ignore();
408  // }
409 
410  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
411  secNumberInit = fvect->size();
412  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
413  secNumberFinal = fvect->size();
414  }
415 
416  G4double secondaryKinetic=-1000*eV;
417 
418  if (!fasterCode) secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
419 
420  // SI - 01/04/2014
421  if (fasterCode)
422  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
423  //
424 
425  G4ThreeVector deltaDirection =
426  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
427  Z, ionizationShell,
428  couple->GetMaterial());
429 
430  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
431 
432  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
433  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
434  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
435  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
436  finalPx /= finalMomentum;
437  finalPy /= finalMomentum;
438  finalPz /= finalMomentum;
439 
440  G4ThreeVector direction;
441  direction.set(finalPx,finalPy,finalPz);
442 
444 
445  // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
446  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
447  G4double deexSecEnergy = 0;
448  for (G4int j=secNumberInit; j < secNumberFinal; j++)
449  {
450  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
451  }
452 
453  if (!statCode)
454  {
456  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
457  }
458  else
459  {
462  }
463 
464  // SI - 01/04/2014
465  if (secondaryKinetic>0)
466  {
467  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
468  fvect->push_back(dp);
469  }
470  //
471 
472  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
474  ionizationShell,
475  theIncomingTrack);
476  }
477 
478 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
double x() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double electron_mass_c2
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAChemistryManager * Instance()
G4double GetPDGMass() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
Hep3Vector unit() const
double y() const
const G4Track * GetCurrentTrack() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4double bindingEnergy(G4int A, G4int Z)
G4AtomicShellEnumerator
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4DNAEmfietzoglouIonisationModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 183 of file G4DNAEmfietzoglouIonisationModel.hh.

184 {
185  fasterCode = input;
186 }
void G4DNAEmfietzoglouIonisationModel::SelectStationary ( G4bool  input)
inline

Definition at line 190 of file G4DNAEmfietzoglouIonisationModel.hh.

191 {
192  statCode = input;
193 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNAEmfietzoglouIonisationModel::fParticleChangeForGamma
protected

Definition at line 91 of file G4DNAEmfietzoglouIonisationModel.hh.


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