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

#include <G4DNABornIonisationModel.hh>

+ Inheritance diagram for G4DNABornIonisationModel:

Public Member Functions

 G4DNABornIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel")
 
virtual ~G4DNABornIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

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

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

Constructor & Destructor Documentation

◆ G4DNABornIonisationModel()

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

Definition at line 43 of file G4DNABornIonisationModel.cc.

45 :G4VEmModel(nam),isInitialised(false)
46{
47 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
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
64 fAtomDeexcitation = 0;
66 fpMolWaterDensity = 0;
67}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4DNABornIonisationModel()

G4DNABornIonisationModel::~G4DNABornIonisationModel ( )
virtual

Definition at line 71 of file G4DNABornIonisationModel.cc.

72{
73 // Cross section
74
75 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
76 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
77 {
78 G4DNACrossSectionDataSet* table = pos->second;
79 delete table;
80 }
81
82 // Final state
83
84 eVecm.clear();
85 pVecm.clear();
86
87}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 243 of file G4DNABornIonisationModel.cc.

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

◆ DifferentialCrossSection()

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

Definition at line 584 of file G4DNABornIonisationModel.cc.

588{
589 G4double sigma = 0.;
590
591 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
592 {
593 G4double valueT1 = 0;
594 G4double valueT2 = 0;
595 G4double valueE21 = 0;
596 G4double valueE22 = 0;
597 G4double valueE12 = 0;
598 G4double valueE11 = 0;
599
600 G4double xs11 = 0;
601 G4double xs12 = 0;
602 G4double xs21 = 0;
603 G4double xs22 = 0;
604
605 if (particleDefinition == G4Electron::ElectronDefinition())
606 {
607 // k should be in eV and energy transfer eV also
608
609 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
610
611 std::vector<double>::iterator t1 = t2-1;
612
613 // SI : the following condition avoids situations where energyTransfer >last vector element
614 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
615 {
616 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
617 std::vector<double>::iterator e11 = e12-1;
618
619 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
620 std::vector<double>::iterator e21 = e22-1;
621
622 valueT1 =*t1;
623 valueT2 =*t2;
624 valueE21 =*e21;
625 valueE22 =*e22;
626 valueE12 =*e12;
627 valueE11 =*e11;
628
629 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
630 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
631 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
632 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
633 }
634
635 }
636
637 if (particleDefinition == G4Proton::ProtonDefinition())
638 {
639 // k should be in eV and energy transfer eV also
640 std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
641 std::vector<double>::iterator t1 = t2-1;
642
643 std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
644 std::vector<double>::iterator e11 = e12-1;
645
646 std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), 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 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
657 xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
658 xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
659 xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
660
661 }
662
663 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
664 if (xsProduct != 0.)
665 {
666 sigma = QuadInterpolator( valueE11, valueE12,
667 valueE21, valueE22,
668 xs11, xs12,
669 xs21, xs22,
670 valueT1, valueT2,
671 k, energyTransfer);
672 }
673
674 }
675
676 return sigma;
677}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 91 of file G4DNABornIonisationModel.cc.

93{
94
95 if (verboseLevel > 3)
96 G4cout << "Calling G4DNABornIonisationModel::Initialise()" << G4endl;
97
98 // Energy limits
99
100 G4String fileElectron("dna/sigma_ionisation_e_born");
101 G4String fileProton("dna/sigma_ionisation_p_born");
102
105
106 G4String electron;
108
109 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
110
111 char *path = getenv("G4LEDATA");
112
113 // *** ELECTRON
114
115 electron = electronDef->GetParticleName();
116
117 tableFile[electron] = fileElectron;
118
119 lowEnergyLimit[electron] = 11. * eV;
120 highEnergyLimit[electron] = 1. * MeV;
121
122 // Cross section
123
125 tableE->LoadData(fileElectron);
126
127 tableData[electron] = tableE;
128
129 // Final state
130
131 std::ostringstream eFullFileName;
132 eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
133 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
134
135 if (!eDiffCrossSection)
136 {
137 G4Exception("G4DNABornIonisationModel::Initialise","em0003",
138 FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_born.dat");
139 }
140
141 eTdummyVec.push_back(0.);
142 while(!eDiffCrossSection.eof())
143 {
144 double tDummy;
145 double eDummy;
146 eDiffCrossSection>>tDummy>>eDummy;
147 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
148 for (int j=0; j<5; j++)
149 {
150 eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
151
152 // SI - only if eof is not reached !
153 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
154
155 eVecm[tDummy].push_back(eDummy);
156
157 }
158 }
159
160 // *** PROTON
161
162 proton = protonDef->GetParticleName();
163
164 tableFile[proton] = fileProton;
165
166 lowEnergyLimit[proton] = 500. * keV;
167 highEnergyLimit[proton] = 100. * MeV;
168
169 // Cross section
170
172 tableP->LoadData(fileProton);
173
174 tableData[proton] = tableP;
175
176 // Final state
177
178 std::ostringstream pFullFileName;
179 pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
180 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
181
182 if (!pDiffCrossSection)
183 {
184 G4Exception("G4DNABornIonisationModel::Initialise","em0003",
185 FatalException,"Missing data file:/dna/sigmadiff_ionisation_p_born.dat");
186 }
187
188 pTdummyVec.push_back(0.);
189 while(!pDiffCrossSection.eof())
190 {
191 double tDummy;
192 double eDummy;
193 pDiffCrossSection>>tDummy>>eDummy;
194 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
195 for (int j=0; j<5; j++)
196 {
197 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
198
199 // SI - only if eof is not reached !
200 if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
201
202 pVecm[tDummy].push_back(eDummy);
203 }
204 }
205
206 //
207
208 if (particle==electronDef)
209 {
210 SetLowEnergyLimit(lowEnergyLimit[electron]);
211 SetHighEnergyLimit(highEnergyLimit[electron]);
212 }
213
214 if (particle==protonDef)
215 {
216 SetLowEnergyLimit(lowEnergyLimit[proton]);
217 SetHighEnergyLimit(highEnergyLimit[proton]);
218 }
219
220 if( verboseLevel>0 )
221 {
222 G4cout << "Born ionisation model is initialized " << G4endl
223 << "Energy range: "
224 << LowEnergyLimit() / eV << " eV - "
225 << HighEnergyLimit() / keV << " keV for "
226 << particle->GetParticleName()
227 << G4endl;
228 }
229
230 // Initialize water density pointer
232
233 //
234 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
235
236 if (isInitialised) { return; }
238 isInitialised = true;
239}
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 326 of file G4DNABornIonisationModel.cc.

331{
332
333 if (verboseLevel > 3)
334 G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel" << G4endl;
335
336 G4double lowLim = 0;
337 G4double highLim = 0;
338
339 G4double k = particle->GetKineticEnergy();
340
341 const G4String& particleName = particle->GetDefinition()->GetParticleName();
342
343 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
344 pos1 = lowEnergyLimit.find(particleName);
345
346 if (pos1 != lowEnergyLimit.end())
347 {
348 lowLim = pos1->second;
349 }
350
351 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
352 pos2 = highEnergyLimit.find(particleName);
353
354 if (pos2 != highEnergyLimit.end())
355 {
356 highLim = pos2->second;
357 }
358
359 if (k >= lowLim && k < highLim)
360 {
361 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
362 G4double particleMass = particle->GetDefinition()->GetPDGMass();
363 G4double totalEnergy = k + particleMass;
364 G4double pSquare = k * (totalEnergy + particleMass);
365 G4double totalMomentum = std::sqrt(pSquare);
366
367 G4int ionizationShell = RandomSelect(k,particleName);
368
369 // sample deexcitation
370 // here we assume that H_{2}O electronic levels are the same of Oxigen.
371 // this can be considered true with a rough 10% error in energy on K-shell,
372
373 G4int secNumberInit = 0; // need to know at a certain point the enrgy of secondaries
374 G4int secNumberFinal = 0; // So I'll make the diference and then sum the energies
375
377 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
378
379 if(fAtomDeexcitation) {
380 G4int Z = 8;
382
383 if (ionizationShell <5 && ionizationShell >1)
384 {
385 as = G4AtomicShellEnumerator(4-ionizationShell);
386 }
387 else if (ionizationShell <2)
388 {
390 }
391
392 // FOR DEBUG ONLY
393 // if (ionizationShell == 4) {
394 //
395 // G4cout << "Z: " << Z << " as: " << as
396 // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
397 // G4cout << "Press <Enter> key to continue..." << G4endl;
398 // G4cin.ignore();
399 // }
400
401 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
402 secNumberInit = fvect->size();
403 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
404 secNumberFinal = fvect->size();
405 }
406
407 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
408
409 G4double cosTheta = 0.;
410 G4double phi = 0.;
411 RandomizeEjectedElectronDirection(particle->GetDefinition(), k,secondaryKinetic, cosTheta, phi);
412
413 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
414 G4double dirX = sinTheta*std::cos(phi);
415 G4double dirY = sinTheta*std::sin(phi);
416 G4double dirZ = cosTheta;
417 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
418 deltaDirection.rotateUz(primaryDirection);
419
421 {
422 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
423
424 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
425 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
426 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
427 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
428 finalPx /= finalMomentum;
429 finalPy /= finalMomentum;
430 finalPz /= finalMomentum;
431
432 G4ThreeVector direction;
433 direction.set(finalPx,finalPy,finalPz);
434
436 }
437
438 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection) ;
439
440 // note thta secondaryKinetic is the nergy of the delta ray, not of all secondaries.
441 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
442 G4double deexSecEnergy = 0;
443 for (G4int j=secNumberInit; j < secNumberFinal; j++) {
444
445 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
446
447 }
448
450 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
451
452 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
453 fvect->push_back(dp);
454
455
456 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
458 ionizationShell,
459 theIncomingTrack);
460 }
461
462}
G4AtomicShellEnumerator
@ eIonizedMolecule
int G4int
Definition: G4Types.hh:66
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 * Electron()
Definition: G4Electron.cc:94
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)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornIonisationModel::fParticleChangeForGamma
protected

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