Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel2 Class Reference

#include <G4DNABornIonisationModel2.hh>

+ Inheritance 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)
 
G4double 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

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

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 *)
 

Detailed Description

Definition at line 47 of file G4DNABornIonisationModel2.hh.

Constructor & Destructor Documentation

◆ G4DNABornIonisationModel2()

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

Definition at line 45 of file G4DNABornIonisationModel2.cc.

46 :
47G4VEmModel(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
65 fAtomDeexcitation = 0;
67 fpMolWaterDensity = 0;
68 fTableData = 0;
69 fLowEnergyLimit = 0;
70 fHighEnergyLimit = 0;
71 fParticleDef = 0;
72
73 // Define default angular generator
74
76
77 // Selection of computation method
78
79 fasterCode = false;
80
81 // Selection of stationary mode
82
83 statCode = false;
84
85 // Selection of SP scaling
86
87 spScaling = true;
88}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4DNABornIonisationModel2()

G4DNABornIonisationModel2::~G4DNABornIonisationModel2 ( )
virtual

Definition at line 92 of file G4DNABornIonisationModel2.cc.

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

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 263 of file G4DNABornIonisationModel2.cc.

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

◆ DifferentialCrossSection()

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

Definition at line 615 of file G4DNABornIonisationModel2.cc.

619{
620 G4double sigma = 0.;
621
622 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)/eV)
623 {
624 G4double valueT1 = 0;
625 G4double valueT2 = 0;
626 G4double valueE21 = 0;
627 G4double valueE22 = 0;
628 G4double valueE12 = 0;
629 G4double valueE11 = 0;
630
631 G4double xs11 = 0;
632 G4double xs12 = 0;
633 G4double xs21 = 0;
634 G4double xs22 = 0;
635
636 // Protection against out of boundary access - proton case : 100 MeV
637 if (k==fTdummyVec.back()) k=k*(1.-1e-12);
638 //
639
640 // k should be in eV and energy transfer eV also
641
642 std::vector<G4double>::iterator t2 = std::upper_bound(fTdummyVec.begin(),
643 fTdummyVec.end(),
644 k);
645
646 std::vector<G4double>::iterator t1 = t2 - 1;
647
648 // SI : the following condition avoids situations where energyTransfer >last vector element
649
650 if (energyTransfer <= fVecm[(*t1)].back()
651 && energyTransfer <= fVecm[(*t2)].back())
652 {
653 std::vector<G4double>::iterator e12 = std::upper_bound(fVecm[(*t1)].begin(),
654 fVecm[(*t1)].end(),
655 energyTransfer);
656 std::vector<G4double>::iterator e11 = e12 - 1;
657
658 std::vector<G4double>::iterator e22 = std::upper_bound(fVecm[(*t2)].begin(),
659 fVecm[(*t2)].end(),
660 energyTransfer);
661 std::vector<G4double>::iterator e21 = e22 - 1;
662
663 valueT1 = *t1;
664 valueT2 = *t2;
665 valueE21 = *e21;
666 valueE22 = *e22;
667 valueE12 = *e12;
668 valueE11 = *e11;
669
670 xs11 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
671 xs12 = fDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
672 xs21 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
673 xs22 = fDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
674
675 }
676
677 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
678 if (xsProduct != 0.)
679 {
680 sigma = QuadInterpolator(valueE11,
681 valueE12,
682 valueE21,
683 valueE22,
684 xs11,
685 xs12,
686 xs21,
687 xs22,
688 valueT1,
689 valueT2,
690 k,
691 energyTransfer);
692 }
693 }
694
695 return sigma;
696}

◆ GetPartialCrossSection()

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

Reimplemented from G4VEmModel.

Definition at line 790 of file G4DNABornIonisationModel2.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 106 of file G4DNABornIonisationModel2.cc.

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

◆ SampleSecondaries()

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 G4double secondaryKinetic=-1000*eV;
353
354 if (fasterCode == false)
355 {
356 secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
357 }
358 else
359 {
360 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
361 }
362
363 G4int Z = 8;
364
365 G4ThreeVector deltaDirection =
366 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
367 Z, ionizationShell,
368 couple->GetMaterial());
369
370 if (secondaryKinetic>0)
371 {
372 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
373 fvect->push_back(dp);
374 }
375
377 {
378 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
379
380 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
381 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
382 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
383 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
384 finalPx /= finalMomentum;
385 finalPy /= finalMomentum;
386 finalPz /= finalMomentum;
387
388 G4ThreeVector direction;
389 direction.set(finalPx,finalPy,finalPz);
390
392 }
393
395
396 // AM: sample deexcitation
397 // here we assume that H_{2}O electronic levels are the same as Oxygen.
398 // this can be considered true with a rough 10% error in energy on K-shell,
399
400 size_t secNumberInit = 0;
401 size_t secNumberFinal = 0;
402
404 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
405
406 // SI: additional protection if tcs interpolation method is modified
407 if (k<bindingEnergy) return;
408 //
409
410 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
411
412 // SI: only atomic deexcitation from K shell is considered
413 if(fAtomDeexcitation && ionizationShell == 4)
414 {
415 const G4AtomicShell* shell =
416 fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
417 secNumberInit = fvect->size();
418 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
419 secNumberFinal = fvect->size();
420
421 if(secNumberFinal > secNumberInit)
422 {
423 for (size_t i=secNumberInit; i<secNumberFinal; ++i)
424 {
425 //Check if there is enough residual energy
426 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
427 {
428 //Ok, this is a valid secondary: keep it
429 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
430 }
431 else
432 {
433 //Invalid secondary: not enough energy to create it!
434 //Keep its energy in the local deposit
435 delete (*fvect)[i];
436 (*fvect)[i]=0;
437 }
438 }
439 }
440
441 }
442
443 //This should never happen
444 if(bindingEnergy < 0.0)
445 G4Exception("G4DNAEmfietzoglouIonisatioModel1::SampleSecondaries()",
446 "em2050",FatalException,"Negative local energy deposit");
447
448 //bindingEnergy has been decreased
449 //by the amount of energy taken away by deexc. products
450 if (!statCode)
451 {
454 }
455 else
456 {
459 }
460
461 // TEST //////////////////////////
462 // if (secondaryKinetic<0) abort();
463 // if (scatteredEnergy<0) abort();
464 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
465 // if (k-scatteredEnergy<0) abort();
466 /////////////////////////////////
467
468 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
470 ionizationShell,
471 theIncomingTrack);
472 }
473}
G4AtomicShellEnumerator
@ eIonizedMolecule
int G4int
Definition: G4Types.hh:85
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4Material * GetMaterial() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectFasterComputation()

void G4DNABornIonisationModel2::SelectFasterComputation ( G4bool  input)
inline

Definition at line 162 of file G4DNABornIonisationModel2.hh.

163{
164 fasterCode = input;
165}

◆ SelectSPScaling()

void G4DNABornIonisationModel2::SelectSPScaling ( G4bool  input)
inline

Definition at line 176 of file G4DNABornIonisationModel2.hh.

177{
178 spScaling = input;
179}

◆ SelectStationary()

void G4DNABornIonisationModel2::SelectStationary ( G4bool  input)
inline

Definition at line 169 of file G4DNABornIonisationModel2.hh.

170{
171 statCode = input;
172}

◆ TransferedEnergy()

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

Definition at line 866 of file G4DNABornIonisationModel2.cc.

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

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornIonisationModel2::fParticleChangeForGamma
protected

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