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

#include <G4DNAEmfietzoglouIonisationModel.hh>

+ Inheritance diagram for G4DNAEmfietzoglouIonisationModel:

Public Member Functions

 G4DNAEmfietzoglouIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
 
virtual ~G4DNAEmfietzoglouIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
G4double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
void SelectFasterComputation (G4bool input)
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void 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 52 of file G4DNAEmfietzoglouIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAEmfietzoglouIonisationModel()

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

Definition at line 52 of file G4DNAEmfietzoglouIonisationModel.cc.

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

◆ ~G4DNAEmfietzoglouIonisationModel()

G4DNAEmfietzoglouIonisationModel::~G4DNAEmfietzoglouIonisationModel ( )
virtual

Definition at line 92 of file G4DNAEmfietzoglouIonisationModel.cc.

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

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 238 of file G4DNAEmfietzoglouIonisationModel.cc.

244{
245 if(verboseLevel > 3)
246 {
247 G4cout
248 << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouIonisationModel"
249 << G4endl;
250 }
251
252 if (particleDefinition != G4Electron::ElectronDefinition()) return 0; // necessary ??
253
254 // Calculate total cross section for model
255
256 G4double sigma=0;
257
258 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
259
260 const G4String& particleName = particleDefinition->GetParticleName();
261
262 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
263 {
264 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
265 pos = tableData.find(particleName);
266
267 if (pos != tableData.end())
268 {
269 G4DNACrossSectionDataSet* table = pos->second;
270 if (table != 0)
271 {
272 sigma = table->FindValue(ekin);
273 }
274 }
275 else
276 {
277 G4Exception("G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume","em0002",
278 FatalException,"Model not applicable to particle type.");
279 }
280 }
281
282 if (verboseLevel > 2)
283 {
284 G4cout << "__________________________________" << G4endl;
285 G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO START" << G4endl;
286 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
287 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
288 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
289 G4cout << "G4DNAEmfietzoglouIonisationModel - XS INFO END" << G4endl;
290 }
291
292 return sigma*waterDensity;
293}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
size_t GetIndex() const
Definition: G4Material.hh:258
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ DifferentialCrossSection()

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

Definition at line 555 of file G4DNAEmfietzoglouIonisationModel.cc.

559{
560 G4double sigma = 0.;
561
562 if(energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex)/eV)
563 {
564 G4double valueT1 = 0;
565 G4double valueT2 = 0;
566 G4double valueE21 = 0;
567 G4double valueE22 = 0;
568 G4double valueE12 = 0;
569 G4double valueE11 = 0;
570
571 G4double xs11 = 0;
572 G4double xs12 = 0;
573 G4double xs21 = 0;
574 G4double xs22 = 0;
575
576 if(particleDefinition == G4Electron::ElectronDefinition())
577 {
578 // Protection against out of boundary access
579 if (k==eTdummyVec.back()) k=k*(1.-1e-12);
580 //
581
582 // k should be in eV and energy transfer eV also
583
584 std::vector<G4double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
585 eTdummyVec.end(),
586 k);
587
588 std::vector<G4double>::iterator t1 = t2 - 1;
589
590 // SI : the following condition avoids situations where energyTransfer >last vector element
591 // added strict limitations (09/08/2017)
592 if(energyTransfer < eVecm[(*t1)].back() &&
593 energyTransfer < eVecm[(*t2)].back())
594 {
595 std::vector<G4double>::iterator e12 =
596 std::upper_bound(eVecm[(*t1)].begin(),
597 eVecm[(*t1)].end(),
598 energyTransfer);
599 std::vector<G4double>::iterator e11 = e12 - 1;
600
601 std::vector<G4double>::iterator e22 =
602 std::upper_bound(eVecm[(*t2)].begin(),
603 eVecm[(*t2)].end(),
604 energyTransfer);
605 std::vector<G4double>::iterator e21 = e22 - 1;
606
607 valueT1 = *t1;
608 valueT2 = *t2;
609 valueE21 = *e21;
610 valueE22 = *e22;
611 valueE12 = *e12;
612 valueE11 = *e11;
613
614 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
615 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
616 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
617 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
618
619 //G4cout << "-------------------" << G4endl;
620 //G4cout << "ionizationLevelIndex=" << ionizationLevelIndex << G4endl;
621 //G4cout << "valueT1/eV=" << valueT1 << " valueT2/eV=" << valueT2 << G4endl;
622 //G4cout << "valueE11/eV=" << valueE11 << " valueE12/eV=" << valueE12
623 // << " valueE21/eV=" << valueE21 << " valueE22/eV=" << valueE22 << G4endl;
624 //G4cout << "xs11=" << xs11 / ((1.e-22 / 3.343) * m*m) << G4endl;
625 //G4cout << "xs12=" << xs12 / ((1.e-22 / 3.343) * m*m) << G4endl;
626 //G4cout << "xs21=" << xs21 / ((1.e-22 / 3.343) * m*m) << G4endl;
627 //G4cout << "xs22=" << xs22 / ((1.e-22 / 3.343) * m*m) << G4endl;
628 //G4cout << "###################" << G4endl;
629
630 }
631
632 }
633
634 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
635 if(xsProduct != 0.)
636 {
637 sigma = QuadInterpolator(valueE11,
638 valueE12,
639 valueE21,
640 valueE22,
641 xs11,
642 xs12,
643 xs21,
644 xs22,
645 valueT1,
646 valueT2,
647 k,
648 energyTransfer);
649 }
650
651 }
652
653 return sigma;
654}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 111 of file G4DNAEmfietzoglouIonisationModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 297 of file G4DNAEmfietzoglouIonisationModel.cc.

303{
304
305 if(verboseLevel > 3)
306 {
307 G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouIonisationModel"
308 << G4endl;
309 }
310
311 G4double k = particle->GetKineticEnergy();
312
313 const G4String& particleName = particle->GetDefinition()->GetParticleName();
314
315 if (k >= LowEnergyLimit() && k <= HighEnergyLimit())
316 {
317 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
318 G4double particleMass = particle->GetDefinition()->GetPDGMass();
319 G4double totalEnergy = k + particleMass;
320 G4double pSquare = k * (totalEnergy + particleMass);
321 G4double totalMomentum = std::sqrt(pSquare);
322
323 G4int ionizationShell = 0;
324
325 ionizationShell = RandomSelect(k,particleName);
326
328 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
329
330 // SI : additional protection if tcs interpolation method is modified
331 if (k<bindingEnergy) return;
332 //
333
334 G4double secondaryKinetic=-1000*eV;
335
336 if (!fasterCode) secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
337
338 if (fasterCode)
339 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
340
341 // SI - For atom. deexc. tagging - 23/05/2017
342
343 G4int Z = 8;
344
345 G4ThreeVector deltaDirection =
346 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
347 Z, ionizationShell,
348 couple->GetMaterial());
349
350 if (secondaryKinetic>0)
351 {
352 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
353 fvect->push_back(dp);
354 }
355
356 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
357
358 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
359 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
360 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
361 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
362 finalPx /= finalMomentum;
363 finalPy /= finalMomentum;
364 finalPz /= finalMomentum;
365
366 G4ThreeVector direction;
367 direction.set(finalPx,finalPy,finalPz);
368
370
371 // AM: sample deexcitation
372 // here we assume that H_{2}O electronic levels are the same as Oxygen.
373 // this can be considered true with a rough 10% error in energy on K-shell,
374
375 size_t secNumberInit = 0;// need to know at a certain point the energy of secondaries
376 size_t secNumberFinal = 0;// So I'll make the diference and then sum the energies
377
378 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
379
380 // SI: only atomic deexcitation from K shell is considered
381 if(fAtomDeexcitation && ionizationShell == 4)
382 {
383 const G4AtomicShell* shell
384 = fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
385 secNumberInit = fvect->size();
386 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
387 secNumberFinal = fvect->size();
388
389 if(secNumberFinal > secNumberInit) {
390 for (size_t i=secNumberInit; i<secNumberFinal; ++i) {
391 //Check if there is enough residual energy
392 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
393 {
394 //Ok, this is a valid secondary: keep it
395 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
396 }
397 else
398 {
399 //Invalid secondary: not enough energy to create it!
400 //Keep its energy in the local deposit
401 delete (*fvect)[i];
402 (*fvect)[i]=0;
403 }
404 }
405 }
406
407 }
408
409 //This should never happen
410 if(bindingEnergy < 0.0)
411 G4Exception("G4DNAEmfietzoglouIonisatioModel1::SampleSecondaries()",
412 "em2050",FatalException,"Negative local energy deposit");
413
414 //bindingEnergy has been decreased
415 //by the amount of energy taken away by deexc. products
416 if (!statCode)
417 {
420 }
421 else
422 {
425 }
426 // TEST //////////////////////////
427 // if (secondaryKinetic<0) abort();
428 // if (scatteredEnergy<0) abort();
429 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
430 // if (k-scatteredEnergy<0) abort();
431 /////////////////////////////////
432
433 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
435 ionizationShell,
436 theIncomingTrack);
437 }
438
439}
G4AtomicShellEnumerator
@ eIonizedMolecule
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: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 G4DNAEmfietzoglouIonisationModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 178 of file G4DNAEmfietzoglouIonisationModel.hh.

179{
180 fasterCode = input;
181}

Referenced by G4EmDNAPhysics_option5::ConstructProcess().

◆ SelectStationary()

void G4DNAEmfietzoglouIonisationModel::SelectStationary ( G4bool  input)
inline

Definition at line 185 of file G4DNAEmfietzoglouIonisationModel.hh.

186{
187 statCode = input;
188}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAEmfietzoglouIonisationModel::fParticleChangeForGamma
protected

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