Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNABornIonisationModel2 Class Reference

#include <G4DNABornIonisationModel2.hh>

Inheritance diagram for G4DNABornIonisationModel2:
Collaboration diagram for G4DNABornIonisationModel2:

Public Member Functions

 G4DNABornIonisationModel2 (const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel")
 
virtual ~G4DNABornIonisationModel2 ()
 
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)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
G4double TransferedEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
 
void SelectFasterComputation (G4bool input)
 
void SelectStationary (G4bool input)
 
void SelectSPScaling (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 * > *)
 
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 SetFluctuationFlag (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
 
G4bool lossFlucFlag
 

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 48 of file G4DNABornIonisationModel2.hh.

Constructor & Destructor Documentation

G4DNABornIonisationModel2::G4DNABornIonisationModel2 ( const G4ParticleDefinition p = 0,
const G4String nam = "DNABornIonisationModel" 
)

Definition at line 46 of file G4DNABornIonisationModel2.cc.

47  :
48  G4VEmModel(nam), isInitialised(false)
49 {
50  verboseLevel = 0;
51  // Verbosity scale:
52  // 0 = nothing
53  // 1 = warning for energy non-conservation
54  // 2 = details of energy budget
55  // 3 = calculation of cross sections, file openings, sampling of atoms
56  // 4 = entering in methods
57 
58  if (verboseLevel > 0)
59  {
60  G4cout << "Born ionisation model is constructed " << G4endl;
61  }
62 
63  //Mark this model as "applicable" for atomic deexcitation
64  SetDeexcitationFlag(true);
65  fAtomDeexcitation = 0;
67  fpMolWaterDensity = 0;
68  fTableData = 0;
69  fLowEnergyLimit = 0;
70  fHighEnergyLimit = 0;
71  fParticleDef = 0;
72 
73  // define default angular generator
75 
76  // Selection of computation method
77 
78  fasterCode = false;
79 
80  // Selection of stationary mode
81 
82  statCode = false;
83 
84  // Selection of SP scaling
85 
86  spScaling = true;
87 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788

Here is the call graph for this function:

G4DNABornIonisationModel2::~G4DNABornIonisationModel2 ( )
virtual

Definition at line 91 of file G4DNABornIonisationModel2.cc.

92 {
93  // Cross section
94 
95  if (fTableData)
96  delete fTableData;
97 
98  // Final state
99 
100  fVecm.clear();
101 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 260 of file G4DNABornIonisationModel2.cc.

265 {
266  if (verboseLevel > 3)
267  {
268  G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel2"
269  << G4endl;
270  }
271 
272  if (particleDefinition != fParticleDef) return 0;
273 
274  // Calculate total cross section for model
275 
276  G4double sigma=0;
277 
278  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
279 
280  if(waterDensity!= 0.0)
281  {
282  if (ekin >= fLowEnergyLimit && ekin < fHighEnergyLimit)
283  {
284  sigma = fTableData->FindValue(ekin);
285 
286  // ICRU49 electronic SP scaling - ZF, SI
287 
288  if (particleDefinition == G4Proton::ProtonDefinition() && ekin < 70*MeV && spScaling)
289  {
290  G4double A = 1.39241700556072800000E-009 ;
291  G4double B = -8.52610412942622630000E-002 ;
292  sigma = sigma * G4Exp(A*(ekin/eV)+B);
293  }
294  //
295  }
296 
297  if (verboseLevel > 2)
298  {
299  G4cout << "__________________________________" << G4endl;
300  G4cout << "G4DNABornIonisationModel2 - XS INFO START" << G4endl;
301  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
302  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
303  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
304  G4cout << "G4DNABornIonisationModel2 - XS INFO END" << G4endl;
305  }
306  } // if (waterMaterial)
307 
308  return sigma*waterDensity;
309 }
size_t GetIndex() const
Definition: G4Material.hh:262
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
double B(double temperature)
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
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
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 607 of file G4DNABornIonisationModel2.cc.

611 {
612  G4double sigma = 0.;
613 
614  if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
615  {
616  G4double valueT1 = 0;
617  G4double valueT2 = 0;
618  G4double valueE21 = 0;
619  G4double valueE22 = 0;
620  G4double valueE12 = 0;
621  G4double valueE11 = 0;
622 
623  G4double xs11 = 0;
624  G4double xs12 = 0;
625  G4double xs21 = 0;
626  G4double xs22 = 0;
627 
628  // k should be in eV and energy transfer eV also
629 
630  std::vector<double>::iterator t2 = std::upper_bound(fTdummyVec.begin(),
631  fTdummyVec.end(),
632  k);
633 
634  std::vector<double>::iterator t1 = t2 - 1;
635 
636  // SI : the following condition avoids situations where energyTransfer >last vector element
637  if (energyTransfer <= fVecm[(*t1)].back()
638  && energyTransfer <= fVecm[(*t2)].back())
639  {
640  std::vector<double>::iterator e12 = std::upper_bound(fVecm[(*t1)].begin(),
641  fVecm[(*t1)].end(),
642  energyTransfer);
643  std::vector<double>::iterator e11 = e12 - 1;
644 
645  std::vector<double>::iterator e22 = std::upper_bound(fVecm[(*t2)].begin(),
646  fVecm[(*t2)].end(),
647  energyTransfer);
648  std::vector<double>::iterator e21 = e22 - 1;
649 
650  valueT1 = *t1;
651  valueT2 = *t2;
652  valueE21 = *e21;
653  valueE22 = *e22;
654  valueE12 = *e12;
655  valueE11 = *e11;
656 
657  xs11 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
658  xs12 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
659  xs21 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
660  xs22 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
661 
662  }
663 
664  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
665  if (xsProduct != 0.)
666  {
667  sigma = QuadInterpolator(valueE11,
668  valueE12,
669  valueE21,
670  valueE22,
671  xs11,
672  xs12,
673  xs21,
674  xs22,
675  valueT1,
676  valueT2,
677  k,
678  energyTransfer);
679  }
680  }
681 
682  return sigma;
683 }
tuple t1
Definition: plottest35.py:33
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 777 of file G4DNABornIonisationModel2.cc.

781 {
782  return fTableData->GetComponent(level)->FindValue(kineticEnergy);
783 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual const G4VEMDataSet * GetComponent(G4int componentId) const

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 105 of file G4DNABornIonisationModel2.cc.

107 {
108 
109  if (verboseLevel > 3)
110  {
111  G4cout << "Calling G4DNABornIonisationModel2::Initialise()" << G4endl;
112  }
113 
114  if(fParticleDef != 0 && particle != fParticleDef)
115  {
116  G4ExceptionDescription description;
117  description << "You are trying to initialized G4DNABornIonisationModel2 "
118  "for particle "
119  << particle->GetParticleName()
120  << G4endl;
121  description << "G4DNABornIonisationModel2 was already initiliased "
122  "for particle:" << fParticleDef->GetParticleName() << G4endl;
123  G4Exception("G4DNABornIonisationModel2::Initialise","bornIonInit",
124  FatalException,description);
125  }
126 
127  fParticleDef = particle;
128 
129  // Energy limits
130  char *path = getenv("G4LEDATA");
131 
132  // ***
133 
134  G4String particleName = particle->GetParticleName();
135  std::ostringstream fullFileName;
136  fullFileName << path;
137 
138  if(particleName == "e-")
139  {
140  fTableFile = "dna/sigma_ionisation_e_born";
141  fLowEnergyLimit = 11.*eV;
142  fHighEnergyLimit = 1.*MeV;
143 
144  if (fasterCode)
145  {
146  fullFileName << "/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat";
147  }
148  else
149  {
150  fullFileName << "/dna/sigmadiff_ionisation_e_born.dat";
151  }
152  }
153  else if(particleName == "proton")
154  {
155  fTableFile = "dna/sigma_ionisation_p_born";
156  fLowEnergyLimit = 500. * keV;
157  fHighEnergyLimit = 100. * MeV;
158 
159  if (fasterCode)
160  {
161  fullFileName << "/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat";
162  }
163  else
164  {
165  fullFileName << "/dna/sigmadiff_ionisation_p_born.dat";
166  }
167  }
168 
169  // Cross section
170  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
171  fTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
172  fTableData->LoadData(fTableFile);
173 
174  // Final state
175 
176  std::ifstream diffCrossSection(fullFileName.str().c_str());
177 
178  if (!diffCrossSection)
179  {
180  G4ExceptionDescription description;
181  description << "Missing data file:" << G4endl << fullFileName.str() << G4endl;
182  G4Exception("G4DNABornIonisationModel2::Initialise","em0003",
183  FatalException,description);
184  }
185 
186  //
187 
188  // Clear the arrays for re-initialization case (MT mode)
189  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
190 
191  fTdummyVec.clear();
192  fVecm.clear();
193 
194  for (int j=0; j<5; j++)
195  {
196  fProbaShellMap[j].clear();
197  fDiffCrossSectionData[j].clear();
198  fNrjTransfData[j].clear();
199  }
200  //
201 
202  fTdummyVec.push_back(0.);
203  while(!diffCrossSection.eof())
204  {
205  double tDummy;
206  double eDummy;
207  diffCrossSection>>tDummy>>eDummy;
208  if (tDummy != fTdummyVec.back()) fTdummyVec.push_back(tDummy);
209 
210  double tmp;
211  for (int j=0; j<5; j++)
212  {
213  diffCrossSection>> tmp;
214 
215  fDiffCrossSectionData[j][tDummy][eDummy] = tmp;
216 
217  if (fasterCode)
218  {
219  fNrjTransfData[j][tDummy][fDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
220  fProbaShellMap[j][tDummy].push_back(fDiffCrossSectionData[j][tDummy][eDummy]);
221  }
222 
223  // SI - only if eof is not reached
224  if (!diffCrossSection.eof() && !fasterCode) fDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
225 
226  if (!fasterCode) fVecm[tDummy].push_back(eDummy);
227 
228  }
229  }
230 
231  //
232  SetLowEnergyLimit(fLowEnergyLimit);
233  SetHighEnergyLimit(fHighEnergyLimit);
234 
235  if( verboseLevel>0 )
236  {
237  G4cout << "Born ionisation model is initialized " << G4endl
238  << "Energy range: "
239  << LowEnergyLimit() / eV << " eV - "
240  << HighEnergyLimit() / keV << " keV for "
241  << particle->GetParticleName()
242  << G4endl;
243  }
244 
245  // Initialize water density pointer
246  fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
247  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
248 
249  //
250  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
251 
252  if (isInitialised)
253  { return;}
255  isInitialised = true;
256 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
static G4LossTableManager * Instance()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
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:732
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double eV
Definition: G4SIunits.hh:215
G4ParticleChangeForGamma * fParticleChangeForGamma
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
static constexpr double MeV
Definition: G4SIunits.hh:214
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 313 of file G4DNABornIonisationModel2.cc.

318 {
319 
320  if (verboseLevel > 3)
321  {
322  G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel2"
323  << G4endl;
324  }
325 
326  G4double k = particle->GetKineticEnergy();
327 
328  if (k >= fLowEnergyLimit && k < fHighEnergyLimit)
329  {
330  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
331  G4double particleMass = particle->GetDefinition()->GetPDGMass();
332  G4double totalEnergy = k + particleMass;
333  G4double pSquare = k * (totalEnergy + particleMass);
334  G4double totalMomentum = std::sqrt(pSquare);
335 
336  G4int ionizationShell = 0;
337 
338  if (!fasterCode) ionizationShell = RandomSelect(k);
339 
340  // SI: The following protection is necessary to avoid infinite loops :
341  // sigmadiff_ionisation_e_born.dat has non zero partial xs at 18 eV for shell 3 (ionizationShell ==2)
342  // sigmadiff_cumulated_ionisation_e_born.dat has zero cumulated partial xs at 18 eV for shell 3 (ionizationShell ==2)
343  // this is due to the fact that the max allowed transfered energy is (18+10.79)/2=17.025 eV and only transfered energies
344  // strictly above this value have non zero partial xs in sigmadiff_ionisation_e_born.dat (starting at trans = 17.12 eV)
345 
346  if (fasterCode)
347  do
348  {
349  ionizationShell = RandomSelect(k);
350  }while (k<19*eV && ionizationShell==2 && particle->GetDefinition()==G4Electron::ElectronDefinition());
351 
352  // AM: sample deexcitation
353  // here we assume that H_{2}O electronic levels are the same as Oxygen.
354  // this can be considered true with a rough 10% error in energy on K-shell,
355 
356  G4int secNumberInit = 0;// need to know at a certain point the energy of secondaries
357  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
358 
360  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
361 
362  //SI: additional protection if tcs interpolation method is modified
363  if (k<bindingEnergy) return;
364  //
365 
366  G4int Z = 8;
367  if(fAtomDeexcitation)
368  {
370 
371  if (ionizationShell <5 && ionizationShell >1)
372  {
373  as = G4AtomicShellEnumerator(4-ionizationShell);
374  }
375  else if (ionizationShell <2)
376  {
377  as = G4AtomicShellEnumerator(3);
378  }
379 
380  // FOR DEBUG ONLY
381  // if (ionizationShell == 4) {
382  //
383  // G4cout << "Z: " << Z << " as: " << as
384  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
385  // G4cout << "Press <Enter> key to continue..." << G4endl;
386  // G4cin.ignore();
387  // }
388 
389  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
390  secNumberInit = fvect->size();
391  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
392  secNumberFinal = fvect->size();
393  }
394 
395  G4double secondaryKinetic=-1000*eV;
396 
397  if (fasterCode == false)
398  {
399  secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
400  }
401  // SI - 01/04/2014
402  else
403  {
404  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
405  }
406  //
407 
408  G4ThreeVector deltaDirection =
409  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
410  Z, ionizationShell,
411  couple->GetMaterial());
412 
413  if (particle->GetDefinition() == G4Electron::ElectronDefinition())
414  {
415  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
416 
417  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
418  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
419  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
420  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
421  finalPx /= finalMomentum;
422  finalPy /= finalMomentum;
423  finalPz /= finalMomentum;
424 
425  G4ThreeVector direction;
426  direction.set(finalPx,finalPy,finalPz);
427 
429  }
430 
431  else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
432 
433  // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
434  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
435  G4double deexSecEnergy = 0;
436  for (G4int j=secNumberInit; j < secNumberFinal; j++)
437  {
438  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
439  }
440 
441  if (!statCode)
442  {
444  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
445  }
446  else
447  {
450  }
451 
452  // SI - 01/04/2014
453  if (secondaryKinetic>0)
454  {
455  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
456  fvect->push_back(dp);
457  }
458  //
459 
460  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
462  ionizationShell,
463  theIncomingTrack);
464  }
465 }
void set(double x, double y, double z)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double GetKineticEnergy() const
double x() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
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
float electron_mass_c2
Definition: hepunit.py:274
G4ParticleChangeForGamma * fParticleChangeForGamma
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 G4DNABornIonisationModel2::SelectFasterComputation ( G4bool  input)
inline

Definition at line 163 of file G4DNABornIonisationModel2.hh.

164 {
165  fasterCode = input;
166 }
void G4DNABornIonisationModel2::SelectSPScaling ( G4bool  input)
inline

Definition at line 177 of file G4DNABornIonisationModel2.hh.

178 {
179  spScaling = input;
180 }
void G4DNABornIonisationModel2::SelectStationary ( G4bool  input)
inline

Definition at line 170 of file G4DNABornIonisationModel2.hh.

171 {
172  statCode = input;
173 }
G4double G4DNABornIonisationModel2::TransferedEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell,
G4double  random 
)

Definition at line 854 of file G4DNABornIonisationModel2.cc.

858 {
859 
860  G4double nrj = 0.;
861 
862  G4double valueK1 = 0;
863  G4double valueK2 = 0;
864  G4double valuePROB21 = 0;
865  G4double valuePROB22 = 0;
866  G4double valuePROB12 = 0;
867  G4double valuePROB11 = 0;
868 
869  G4double nrjTransf11 = 0;
870  G4double nrjTransf12 = 0;
871  G4double nrjTransf21 = 0;
872  G4double nrjTransf22 = 0;
873 
874  // k should be in eV
875  std::vector<double>::iterator k2 = std::upper_bound(fTdummyVec.begin(),
876  fTdummyVec.end(),
877  k);
878  std::vector<double>::iterator k1 = k2 - 1;
879 
880  /*
881  G4cout << "----> k=" << k
882  << " " << *k1
883  << " " << *k2
884  << " " << random
885  << " " << ionizationLevelIndex
886  << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
887  << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
888  << G4endl;
889  */
890 
891  // SI : the following condition avoids situations where random >last vector element
892  if (random <= fProbaShellMap[ionizationLevelIndex][(*k1)].back()
893  && random <= fProbaShellMap[ionizationLevelIndex][(*k2)].back())
894  {
895  std::vector<double>::iterator prob12 =
896  std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
897  fProbaShellMap[ionizationLevelIndex][(*k1)].end(),
898  random);
899 
900  std::vector<double>::iterator prob11 = prob12 - 1;
901 
902  std::vector<double>::iterator prob22 =
903  std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
904  fProbaShellMap[ionizationLevelIndex][(*k2)].end(),
905  random);
906 
907  std::vector<double>::iterator prob21 = prob22 - 1;
908 
909  valueK1 = *k1;
910  valueK2 = *k2;
911  valuePROB21 = *prob21;
912  valuePROB22 = *prob22;
913  valuePROB12 = *prob12;
914  valuePROB11 = *prob11;
915 
916  /*
917  G4cout << " " << random << " " << valuePROB11 << " "
918  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
919  */
920 
921  nrjTransf11 = fNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
922  nrjTransf12 = fNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
923  nrjTransf21 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
924  nrjTransf22 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
925 
926  /*
927  G4cout << " " << ionizationLevelIndex << " "
928  << random << " " <<valueK1 << " " << valueK2 << G4endl;
929 
930  G4cout << " " << random << " " << nrjTransf11 << " "
931  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
932  */
933 
934  }
935  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
936  if (random > fProbaShellMap[ionizationLevelIndex][(*k1)].back())
937  {
938  std::vector<double>::iterator prob22 =
939  std::upper_bound(fProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
940  fProbaShellMap[ionizationLevelIndex][(*k2)].end(),
941  random);
942 
943  std::vector<double>::iterator prob21 = prob22 - 1;
944 
945  valueK1 = *k1;
946  valueK2 = *k2;
947  valuePROB21 = *prob21;
948  valuePROB22 = *prob22;
949 
950  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
951 
952  nrjTransf21 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
953  nrjTransf22 = fNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
954 
955  G4double interpolatedvalue2 = Interpolate(valuePROB21,
956  valuePROB22,
957  random,
958  nrjTransf21,
959  nrjTransf22);
960 
961  // zeros are explicitly set
962 
963  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
964 
965  /*
966  G4cout << " " << ionizationLevelIndex << " "
967  << random << " " <<valueK1 << " " << valueK2 << G4endl;
968 
969  G4cout << " " << random << " " << nrjTransf11 << " "
970  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
971 
972  G4cout << "ici" << " " << value << G4endl;
973  */
974 
975  return value;
976  }
977 
978  // End electron and proton cases
979 
980  G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21
981  * nrjTransf22;
982  //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
983 
984  if (nrjTransfProduct != 0.)
985  {
986  nrj = QuadInterpolator(valuePROB11,
987  valuePROB12,
988  valuePROB21,
989  valuePROB22,
990  nrjTransf11,
991  nrjTransf12,
992  nrjTransf21,
993  nrjTransf22,
994  valueK1,
995  valueK2,
996  k,
997  random);
998  }
999  //G4cout << nrj << endl;
1000 
1001  return nrj;
1002 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
double G4double
Definition: G4Types.hh:76

Member Data Documentation

G4ParticleChangeForGamma* G4DNABornIonisationModel2::fParticleChangeForGamma
protected

Definition at line 89 of file G4DNABornIonisationModel2.hh.


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