Geant4  10.02.p03
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 *> *)
 
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 Types

typedef std::map< double, std::map< double, double > > TriDimensionMap
 
typedef std::map< double, std::vector< double > > VecMap
 

Private Member Functions

G4double RandomizeEjectedElectronEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
 
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
 
G4double Interpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 
G4double QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
 
G4int RandomSelect (G4double energy)
 
G4DNABornIonisationModel2operator= (const G4DNABornIonisationModel2 &right)
 
 G4DNABornIonisationModel2 (const G4DNABornIonisationModel2 &)
 

Private Attributes

G4bool fasterCode
 
G4bool statCode
 
G4bool spScaling
 
const std::vector< G4double > * fpMolWaterDensity
 
G4VAtomDeexcitationfAtomDeexcitation
 
G4double fLowEnergyLimit
 
G4double fHighEnergyLimit
 
const G4ParticleDefinitionfParticleDef
 
G4bool isInitialised
 
G4int verboseLevel
 
G4String fTableFile
 
G4DNACrossSectionDataSetfTableData
 
G4DNAWaterIonisationStructure waterStructure
 
TriDimensionMap fDiffCrossSectionData [6]
 
TriDimensionMap fNrjTransfData [6]
 
std::vector< double > fTdummyVec
 
VecMap fVecm
 
VecMap fProbaShellMap [6]
 

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

Member Typedef Documentation

◆ TriDimensionMap

typedef std::map<double, std::map<double, double> > G4DNABornIonisationModel2::TriDimensionMap
private

Definition at line 138 of file G4DNABornIonisationModel2.hh.

◆ VecMap

typedef std::map<double, std::vector<double> > G4DNABornIonisationModel2::VecMap
private

Definition at line 145 of file G4DNABornIonisationModel2.hh.

Constructor & Destructor Documentation

◆ G4DNABornIonisationModel2() [1/2]

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

Definition at line 45 of file G4DNABornIonisationModel2.cc.

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

◆ ~G4DNABornIonisationModel2()

G4DNABornIonisationModel2::~G4DNABornIonisationModel2 ( )
virtual

Definition at line 90 of file G4DNABornIonisationModel2.cc.

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

◆ G4DNABornIonisationModel2() [2/2]

G4DNABornIonisationModel2::G4DNABornIonisationModel2 ( const G4DNABornIonisationModel2 )
private

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 259 of file G4DNABornIonisationModel2.cc.

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

◆ DifferentialCrossSection()

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

Definition at line 606 of file G4DNABornIonisationModel2.cc.

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

◆ GetPartialCrossSection()

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

Reimplemented from G4VEmModel.

Definition at line 776 of file G4DNABornIonisationModel2.cc.

780 {
781  return fTableData->GetComponent(level)->FindValue(kineticEnergy);
782 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
G4DNACrossSectionDataSet * fTableData
Here is the call graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 104 of file G4DNABornIonisationModel2.cc.

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

◆ Interpolate()

G4double G4DNABornIonisationModel2::Interpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 686 of file G4DNABornIonisationModel2.cc.

691 {
692  G4double value = 0.;
693 
694  // Log-log interpolation by default
695 
696  if (e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
697  && !fasterCode)
698  {
699  G4double a = (std::log10(xs2) - std::log10(xs1))
700  / (std::log10(e2) - std::log10(e1));
701  G4double b = std::log10(xs2) - a * std::log10(e2);
702  G4double sigma = a * std::log10(e) + b;
703  value = (std::pow(10., sigma));
704  }
705 
706  // Switch to lin-lin interpolation
707  /*
708  if ((e2-e1)!=0)
709  {
710  G4double d1 = xs1;
711  G4double d2 = xs2;
712  value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
713  }
714  */
715 
716  // Switch to log-lin interpolation for faster code
717  if ((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
718  {
719  G4double d1 = std::log10(xs1);
720  G4double d2 = std::log10(xs2);
721  value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
722  }
723 
724  // Switch to lin-lin interpolation for faster code
725  // in case one of xs1 or xs2 (=cum proba) value is zero
726 
727  if ((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
728  {
729  G4double d1 = xs1;
730  G4double d2 = xs2;
731  value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
732  }
733 
734  /*
735  G4cout
736  << e1 << " "
737  << e2 << " "
738  << e << " "
739  << xs1 << " "
740  << xs2 << " "
741  << value
742  << G4endl;
743  */
744 
745  return value;
746 }
static const G4double d1
static const G4double e2
static const G4double e1
double G4double
Definition: G4Types.hh:76
static const G4double d2
Here is the caller graph for this function:

◆ operator=()

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

◆ QuadInterpolator()

G4double G4DNABornIonisationModel2::QuadInterpolator ( G4double  e11,
G4double  e12,
G4double  e21,
G4double  e22,
G4double  x11,
G4double  x12,
G4double  x21,
G4double  x22,
G4double  t1,
G4double  t2,
G4double  t,
G4double  e 
)
private

Definition at line 750 of file G4DNABornIonisationModel2.cc.

762 {
763  G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
764  G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
765  G4double value = Interpolate(t1,
766  t2,
767  t,
768  interpolatedvalue1,
769  interpolatedvalue2);
770 
771  return value;
772 }
TTree * t1
Definition: plottest35.C:26
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomizeEjectedElectronEnergy()

G4double G4DNABornIonisationModel2::RandomizeEjectedElectronEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell 
)
private

Definition at line 468 of file G4DNABornIonisationModel2.cc.

471 {
472  // G4cout << "*** SLOW computation for " << " " << particleDefinition->GetParticleName() << G4endl;
473 
474  if (particleDefinition == G4Electron::ElectronDefinition())
475  {
476  G4double maximumEnergyTransfer = 0.;
477  if ((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
478  maximumEnergyTransfer = k;
479  else
480  maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell)) / 2.;
481 
482  // SI : original method
483  /*
484  G4double crossSectionMaximum = 0.;
485  for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
486  {
487  G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
488  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
489  }
490  */
491 
492  // SI : alternative method
493  G4double crossSectionMaximum = 0.;
494 
495  G4double minEnergy = waterStructure.IonisationEnergy(shell);
496  G4double maxEnergy = maximumEnergyTransfer;
497  G4int nEnergySteps = 50;
498 
499  G4double value(minEnergy);
500  G4double stpEnergy(std::pow(maxEnergy / value,
501  1. / static_cast<G4double>(nEnergySteps - 1)));
502  G4int step(nEnergySteps);
503  while (step > 0)
504  {
505  step--;
506  G4double differentialCrossSection =
507  DifferentialCrossSection(particleDefinition,
508  k / eV,
509  value / eV,
510  shell);
511  if (differentialCrossSection >= crossSectionMaximum)
512  crossSectionMaximum = differentialCrossSection;
513  value *= stpEnergy;
514  }
515  //
516 
517  G4double secondaryElectronKineticEnergy = 0.;
518  do
519  {
520  secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
521  }while(G4UniformRand()*crossSectionMaximum >
522  DifferentialCrossSection(particleDefinition, k/eV,
523  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
524 
525  return secondaryElectronKineticEnergy;
526 
527  }
528 
529  else if (particleDefinition == G4Proton::ProtonDefinition())
530  {
531  G4double maximumKineticEnergyTransfer = 4.
533 
534  G4double crossSectionMaximum = 0.;
535  for (G4double value = waterStructure.IonisationEnergy(shell);
536  value <= 4. * waterStructure.IonisationEnergy(shell); value += 0.1 * eV)
537  {
538  G4double differentialCrossSection =
539  DifferentialCrossSection(particleDefinition,
540  k / eV,
541  value / eV,
542  shell);
543  if (differentialCrossSection >= crossSectionMaximum)
544  crossSectionMaximum = differentialCrossSection;
545  }
546 
547  G4double secondaryElectronKineticEnergy = 0.;
548  do
549  {
550  secondaryElectronKineticEnergy = G4UniformRand()* maximumKineticEnergyTransfer;
551  }while(G4UniformRand()*crossSectionMaximum >=
552  DifferentialCrossSection(particleDefinition, k/eV,
553  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
554 
555  return secondaryElectronKineticEnergy;
556  }
557 
558  return 0;
559 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
float proton_mass_c2
Definition: hepunit.py:275
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
float electron_mass_c2
Definition: hepunit.py:274
static const double eV
Definition: G4SIunits.hh:212
G4DNAWaterIonisationStructure waterStructure
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomizeEjectedElectronEnergyFromCumulatedDcs()

G4double G4DNABornIonisationModel2::RandomizeEjectedElectronEnergyFromCumulatedDcs ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell 
)
private

Definition at line 826 of file G4DNABornIonisationModel2.cc.

829 {
830  //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
831 
832  G4double secondaryElectronKineticEnergy = 0.;
833 
834  G4double random = G4UniformRand();
835 
836  secondaryElectronKineticEnergy = TransferedEnergy(particleDefinition,
837  k / eV,
838  shell,
839  random) * eV
840  - waterStructure.IonisationEnergy(shell);
841 
842  //G4cout << TransferedEnergy(particleDefinition, k/eV, shell, random) << G4endl;
843  // SI - 01/04/2014
844  if (secondaryElectronKineticEnergy < 0.)
845  return 0.;
846  //
847 
848  return secondaryElectronKineticEnergy;
849 }
#define G4UniformRand()
Definition: Randomize.hh:97
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
static const double eV
Definition: G4SIunits.hh:212
G4DNAWaterIonisationStructure waterStructure
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomSelect()

G4int G4DNABornIonisationModel2::RandomSelect ( G4double  energy)
private

Definition at line 786 of file G4DNABornIonisationModel2.cc.

787 {
788  G4int level = 0;
789 
790  G4double* valuesBuffer = new G4double[fTableData->NumberOfComponents()];
791  const size_t n(fTableData->NumberOfComponents());
792  size_t i(n);
793  G4double value = 0.;
794 
795  while (i > 0)
796  {
797  i--;
798  valuesBuffer[i] = fTableData->GetComponent(i)->FindValue(k);
799  value += valuesBuffer[i];
800  }
801 
802  value *= G4UniformRand();
803 
804  i = n;
805 
806  while (i > 0)
807  {
808  i--;
809 
810  if (valuesBuffer[i] > value)
811  {
812  delete[] valuesBuffer;
813  return i;
814  }
815  value -= valuesBuffer[i];
816  }
817 
818  if (valuesBuffer)
819  delete[] valuesBuffer;
820 
821  return level;
822 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
int G4int
Definition: G4Types.hh:78
Char_t n[5]
#define G4UniformRand()
Definition: Randomize.hh:97
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
G4DNACrossSectionDataSet * fTableData
double G4double
Definition: G4Types.hh:76
virtual size_t NumberOfComponents(void) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 312 of file G4DNABornIonisationModel2.cc.

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

◆ SelectFasterComputation()

void G4DNABornIonisationModel2::SelectFasterComputation ( G4bool  input)
inline

Definition at line 163 of file G4DNABornIonisationModel2.hh.

164 {
165  fasterCode = input;
166 }

◆ SelectSPScaling()

void G4DNABornIonisationModel2::SelectSPScaling ( G4bool  input)
inline

Definition at line 177 of file G4DNABornIonisationModel2.hh.

178 {
179  spScaling = input;
180 }

◆ SelectStationary()

void G4DNABornIonisationModel2::SelectStationary ( G4bool  input)
inline

Definition at line 170 of file G4DNABornIonisationModel2.hh.

171 {
172  statCode = input;
173 }

◆ TransferedEnergy()

G4double G4DNABornIonisationModel2::TransferedEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell,
G4double  random 
)

Definition at line 853 of file G4DNABornIonisationModel2.cc.

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

Member Data Documentation

◆ fasterCode

G4bool G4DNABornIonisationModel2::fasterCode
private

Definition at line 93 of file G4DNABornIonisationModel2.hh.

◆ fAtomDeexcitation

G4VAtomDeexcitation* G4DNABornIonisationModel2::fAtomDeexcitation
private

Definition at line 101 of file G4DNABornIonisationModel2.hh.

◆ fDiffCrossSectionData

TriDimensionMap G4DNABornIonisationModel2::fDiffCrossSectionData[6]
private

Definition at line 140 of file G4DNABornIonisationModel2.hh.

◆ fHighEnergyLimit

G4double G4DNABornIonisationModel2::fHighEnergyLimit
private

Definition at line 104 of file G4DNABornIonisationModel2.hh.

◆ fLowEnergyLimit

G4double G4DNABornIonisationModel2::fLowEnergyLimit
private

Definition at line 103 of file G4DNABornIonisationModel2.hh.

◆ fNrjTransfData

TriDimensionMap G4DNABornIonisationModel2::fNrjTransfData[6]
private

Definition at line 141 of file G4DNABornIonisationModel2.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornIonisationModel2::fParticleChangeForGamma
protected

Definition at line 89 of file G4DNABornIonisationModel2.hh.

◆ fParticleDef

const G4ParticleDefinition* G4DNABornIonisationModel2::fParticleDef
private

Definition at line 106 of file G4DNABornIonisationModel2.hh.

◆ fpMolWaterDensity

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

Definition at line 98 of file G4DNABornIonisationModel2.hh.

◆ fProbaShellMap

VecMap G4DNABornIonisationModel2::fProbaShellMap[6]
private

Definition at line 148 of file G4DNABornIonisationModel2.hh.

◆ fTableData

G4DNACrossSectionDataSet* G4DNABornIonisationModel2::fTableData
private

Definition at line 113 of file G4DNABornIonisationModel2.hh.

◆ fTableFile

G4String G4DNABornIonisationModel2::fTableFile
private

Definition at line 112 of file G4DNABornIonisationModel2.hh.

◆ fTdummyVec

std::vector<double> G4DNABornIonisationModel2::fTdummyVec
private

Definition at line 143 of file G4DNABornIonisationModel2.hh.

◆ fVecm

VecMap G4DNABornIonisationModel2::fVecm
private

Definition at line 147 of file G4DNABornIonisationModel2.hh.

◆ isInitialised

G4bool G4DNABornIonisationModel2::isInitialised
private

Definition at line 108 of file G4DNABornIonisationModel2.hh.

◆ spScaling

G4bool G4DNABornIonisationModel2::spScaling
private

Definition at line 95 of file G4DNABornIonisationModel2.hh.

◆ statCode

G4bool G4DNABornIonisationModel2::statCode
private

Definition at line 94 of file G4DNABornIonisationModel2.hh.

◆ verboseLevel

G4int G4DNABornIonisationModel2::verboseLevel
private

Definition at line 109 of file G4DNABornIonisationModel2.hh.

◆ waterStructure

G4DNAWaterIonisationStructure G4DNABornIonisationModel2::waterStructure
private

Definition at line 117 of file G4DNABornIonisationModel2.hh.


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