source: trunk/source/processes/electromagnetic/lowenergy/src/G4LowEnergyIonisation.cc@ 830

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 24.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4LowEnergyIonisation.cc,v 1.102 2006/06/29 19:40:19 gunter Exp $
27// GEANT4 tag $Name: $
28//
29// --------------------------------------------------------------
30//
31// File name: G4LowEnergyIonisation
32//
33// Author: Alessandra Forti, Vladimir Ivanchenko
34//
35// Creation date: March 1999
36//
37// Modifications:
38// - 11.04.2000 VL
39// Changing use of float and G4float casts to G4double casts
40// because of problems with optimisation (bug ?)
41// 10.04.2000 VL
42// - Correcting Fluorescence transition probabilities in order to take into account
43// non-radiative transitions. No Auger electron simulated yet: energy is locally deposited.
44// 10.04.2000 VL
45// - Correction of incident electron final momentum direction
46// 07.04.2000 VL+LU
47// - First implementation of continuous energy loss
48// 22.03.2000 VL
49// - 1 bug corrected in SelectRandomAtom method (units)
50// 17.02.2000 Veronique Lefebure
51// - 5 bugs corrected:
52// *in Fluorescence, 2 bugs affecting
53// . localEnergyDeposition and
54// . number of emitted photons that was then always 1 less
55// *in EnergySampling method:
56// . expon = Parms[13]+1; (instead of uncorrect -1)
57// . rejection /= Parms[6];(instead of uncorrect Parms[7])
58// . Parms[6] is apparently corrupted in the data file (often = 0)
59// -->Compute normalisation into local variable rejectionMax
60// and use rejectionMax in stead of Parms[6]
61//
62// Added Livermore data table construction methods A. Forti
63// Modified BuildMeanFreePath to read new data tables A. Forti
64// Added EnergySampling method A. Forti
65// Modified PostStepDoIt to insert sampling with EEDL data A. Forti
66// Added SelectRandomAtom A. Forti
67// Added map of the elements A. Forti
68// 20.09.00 V.Ivanchenko update fluctuations
69// 24.04.01 V.Ivanchenko remove RogueWave
70// 22.05.01 V.Ivanchenko update calculation of delta-ray kinematic +
71// clean up the code
72// 02.08.01 V.Ivanchenko fix energy conservation for small steps
73// 18.08.01 V.Ivanchenko fix energy conservation for pathalogical delta-energy
74// 01.10.01 E. Guardincerri Replaced fluorescence generation in PostStepDoIt
75// according to design iteration
76// 04.10.01 MGP Minor clean-up in the fluo section, removal of
77// compilation warnings and extra protection to
78// prevent from accessing a null pointer
79// 29.09.01 V.Ivanchenko revision based on design iteration
80// 10.10.01 MGP Revision to improve code quality and
81// consistency with design
82// 18.10.01 V.Ivanchenko Add fluorescence AlongStepDoIt
83// 18.10.01 MGP Revision to improve code quality and
84// consistency with design
85// 19.10.01 V.Ivanchenko update according to new design, V.Ivanchenko
86// 26.10.01 V.Ivanchenko clean up deexcitation
87// 28.10.01 V.Ivanchenko update printout
88// 29.11.01 V.Ivanchenko New parametrisation introduced
89// 25.03.02 V.Ivanchneko Fix in fluorescence
90// 28.03.02 V.Ivanchenko Add flag of fluorescence
91// 28.05.02 V.Ivanchenko Remove flag fStopAndKill
92// 31.05.02 V.Ivanchenko Add path of Fluo + Auger cuts to
93// AtomicDeexcitation
94// 03.06.02 MGP Restore fStopAndKill
95// 19.06.02 VI Additional printout
96// 30.07.02 VI Fix in restricted energy loss
97// 20.09.02 VI Remove ActivateFlurescence from SetCut...
98// 21.01.03 VI Cut per region
99// 12.02.03 VI Change signature for Deexcitation
100// 12.04.03 V.Ivanchenko Cut per region for fluo AlongStep
101// 31.08.04 V.Ivanchenko Add density correction
102//
103// --------------------------------------------------------------
104
105#include "G4LowEnergyIonisation.hh"
106#include "G4eIonisationSpectrum.hh"
107#include "G4eIonisationCrossSectionHandler.hh"
108#include "G4AtomicTransitionManager.hh"
109#include "G4AtomicShell.hh"
110#include "G4VDataSetAlgorithm.hh"
111#include "G4SemiLogInterpolation.hh"
112#include "G4LogLogInterpolation.hh"
113#include "G4EMDataSet.hh"
114#include "G4VEMDataSet.hh"
115#include "G4CompositeEMDataSet.hh"
116#include "G4EnergyLossTables.hh"
117#include "G4ShellVacancy.hh"
118#include "G4UnitsTable.hh"
119#include "G4Electron.hh"
120#include "G4Gamma.hh"
121#include "G4ProductionCutsTable.hh"
122
123G4LowEnergyIonisation::G4LowEnergyIonisation(const G4String& nam)
124 : G4eLowEnergyLoss(nam),
125 crossSectionHandler(0),
126 theMeanFreePath(0),
127 energySpectrum(0),
128 shellVacancy(0)
129{
130 cutForPhotons = 250.0*eV;
131 cutForElectrons = 250.0*eV;
132 verboseLevel = 0;
133}
134
135
136G4LowEnergyIonisation::~G4LowEnergyIonisation()
137{
138 delete crossSectionHandler;
139 delete energySpectrum;
140 delete theMeanFreePath;
141 delete shellVacancy;
142}
143
144
145void G4LowEnergyIonisation::BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
146{
147 if(verboseLevel > 0) {
148 G4cout << "G4LowEnergyIonisation::BuildPhysicsTable start"
149 << G4endl;
150 }
151
152 cutForDelta.clear();
153
154 // Create and fill IonisationParameters once
155 if( energySpectrum != 0 ) delete energySpectrum;
156 energySpectrum = new G4eIonisationSpectrum();
157
158 if(verboseLevel > 0) {
159 G4cout << "G4VEnergySpectrum is initialized"
160 << G4endl;
161 }
162
163 // Create and fill G4CrossSectionHandler once
164
165 if ( crossSectionHandler != 0 ) delete crossSectionHandler;
166 G4VDataSetAlgorithm* interpolation = new G4SemiLogInterpolation();
167 G4double lowKineticEnergy = GetLowerBoundEloss();
168 G4double highKineticEnergy = GetUpperBoundEloss();
169 G4int totBin = GetNbinEloss();
170 crossSectionHandler = new G4eIonisationCrossSectionHandler(energySpectrum,
171 interpolation,
172 lowKineticEnergy,
173 highKineticEnergy,
174 totBin);
175 crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
176
177 if (verboseLevel > 0) {
178 G4cout << GetProcessName()
179 << " is created; Cross section data: "
180 << G4endl;
181 crossSectionHandler->PrintData();
182 G4cout << "Parameters: "
183 << G4endl;
184 energySpectrum->PrintData();
185 }
186
187 // Build loss table for IonisationIV
188
189 BuildLossTable(aParticleType);
190
191 if(verboseLevel > 0) {
192 G4cout << "The loss table is built"
193 << G4endl;
194 }
195
196 if (&aParticleType==G4Electron::Electron()) {
197
198 RecorderOfElectronProcess[CounterOfElectronProcess] = (*this).theLossTable;
199 CounterOfElectronProcess++;
200 PrintInfoDefinition();
201
202 } else {
203
204 RecorderOfPositronProcess[CounterOfPositronProcess] = (*this).theLossTable;
205 CounterOfPositronProcess++;
206 }
207
208 // Build mean free path data using cut values
209
210 if( theMeanFreePath ) delete theMeanFreePath;
211 theMeanFreePath = crossSectionHandler->
212 BuildMeanFreePathForMaterials(&cutForDelta);
213
214 if(verboseLevel > 0) {
215 G4cout << "The MeanFreePath table is built"
216 << G4endl;
217 if(verboseLevel > 1) theMeanFreePath->PrintData();
218 }
219
220 // Build common DEDX table for all ionisation processes
221
222 BuildDEDXTable(aParticleType);
223
224 if (verboseLevel > 0) {
225 G4cout << "G4LowEnergyIonisation::BuildPhysicsTable end"
226 << G4endl;
227 }
228}
229
230
231void G4LowEnergyIonisation::BuildLossTable(const G4ParticleDefinition& )
232{
233 // Build table for energy loss due to soft brems
234 // the tables are built for *MATERIALS* binning is taken from LowEnergyLoss
235
236 G4double lowKineticEnergy = GetLowerBoundEloss();
237 G4double highKineticEnergy = GetUpperBoundEloss();
238 size_t totBin = GetNbinEloss();
239
240 // create table
241
242 if (theLossTable) {
243 theLossTable->clearAndDestroy();
244 delete theLossTable;
245 }
246 const G4ProductionCutsTable* theCoupleTable=
247 G4ProductionCutsTable::GetProductionCutsTable();
248 size_t numOfCouples = theCoupleTable->GetTableSize();
249 theLossTable = new G4PhysicsTable(numOfCouples);
250
251 if (shellVacancy != 0) delete shellVacancy;
252 shellVacancy = new G4ShellVacancy();
253 G4DataVector* ksi = 0;
254 G4DataVector* energy = 0;
255 size_t binForFluo = totBin/10;
256
257 G4PhysicsLogVector* bVector = new G4PhysicsLogVector(lowKineticEnergy,
258 highKineticEnergy,
259 binForFluo);
260 const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
261
262 // Clean up the vector of cuts
263
264 cutForDelta.clear();
265
266 // Loop for materials
267
268 for (size_t m=0; m<numOfCouples; m++) {
269
270 // create physics vector and fill it
271 G4PhysicsLogVector* aVector = new G4PhysicsLogVector(lowKineticEnergy,
272 highKineticEnergy,
273 totBin);
274
275 // get material parameters needed for the energy loss calculation
276 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
277 const G4Material* material= couple->GetMaterial();
278
279 // the cut cannot be below lowest limit
280 G4double tCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[m];
281 if(tCut > highKineticEnergy) tCut = highKineticEnergy;
282 cutForDelta.push_back(tCut);
283 const G4ElementVector* theElementVector = material->GetElementVector();
284 size_t NumberOfElements = material->GetNumberOfElements() ;
285 const G4double* theAtomicNumDensityVector =
286 material->GetAtomicNumDensityVector();
287 if(verboseLevel > 0) {
288 G4cout << "Energy loss for material # " << m
289 << " tCut(keV)= " << tCut/keV
290 << G4endl;
291 }
292
293 // now comes the loop for the kinetic energy values
294 for (size_t i = 0; i<totBin; i++) {
295
296 G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(i);
297 G4double ionloss = 0.;
298
299 // loop for elements in the material
300 for (size_t iel=0; iel<NumberOfElements; iel++ ) {
301
302 G4int Z = (G4int)((*theElementVector)[iel]->GetZ());
303
304 G4int nShells = transitionManager->NumberOfShells(Z);
305
306 for (G4int n=0; n<nShells; n++) {
307
308 G4double e = energySpectrum->AverageEnergy(Z, 0.0, tCut,
309 lowEdgeEnergy, n);
310 G4double cs= crossSectionHandler->FindValue(Z, lowEdgeEnergy, n);
311 ionloss += e * cs * theAtomicNumDensityVector[iel];
312
313 if(verboseLevel > 1 || (Z == 14 && lowEdgeEnergy>1. && lowEdgeEnergy<0.)) {
314 G4cout << "Z= " << Z
315 << " shell= " << n
316 << " E(keV)= " << lowEdgeEnergy/keV
317 << " Eav(keV)= " << e/keV
318 << " cs= " << cs
319 << " loss= " << ionloss
320 << " rho= " << theAtomicNumDensityVector[iel]
321 << G4endl;
322 }
323 }
324 G4double esp = energySpectrum->Excitation(Z, lowEdgeEnergy);
325 ionloss += esp * theAtomicNumDensityVector[iel];
326
327 }
328 if(verboseLevel > 1 || (m == 0 && lowEdgeEnergy>=1. && lowEdgeEnergy<=0.)) {
329 G4cout << "Sum: "
330 << " E(keV)= " << lowEdgeEnergy/keV
331 << " loss(MeV/mm)= " << ionloss*mm/MeV
332 << G4endl;
333 }
334 aVector->PutValue(i,ionloss);
335 }
336 theLossTable->insert(aVector);
337
338 // fill data for fluorescence
339
340 G4VDataSetAlgorithm* interp = new G4LogLogInterpolation();
341 G4VEMDataSet* xsis = new G4CompositeEMDataSet(interp, 1., 1.);
342 for (size_t iel=0; iel<NumberOfElements; iel++ ) {
343
344 G4int Z = (G4int)((*theElementVector)[iel]->GetZ());
345 energy = new G4DataVector();
346 ksi = new G4DataVector();
347
348 for (size_t j = 0; j<binForFluo; j++) {
349
350 G4double lowEdgeEnergy = bVector->GetLowEdgeEnergy(j);
351 G4double cross = 0.;
352 G4double eAverage= 0.;
353 G4int nShells = transitionManager->NumberOfShells(Z);
354
355 for (G4int n=0; n<nShells; n++) {
356
357 G4double e = energySpectrum->AverageEnergy(Z, 0.0, tCut,
358 lowEdgeEnergy, n);
359 G4double pro = energySpectrum->Probability(Z, 0.0, tCut,
360 lowEdgeEnergy, n);
361 G4double cs= crossSectionHandler->FindValue(Z, lowEdgeEnergy, n);
362 eAverage += e * cs * theAtomicNumDensityVector[iel];
363 cross += cs * pro * theAtomicNumDensityVector[iel];
364 if(verboseLevel > 1) {
365 G4cout << "Z= " << Z
366 << " shell= " << n
367 << " E(keV)= " << lowEdgeEnergy/keV
368 << " Eav(keV)= " << e/keV
369 << " pro= " << pro
370 << " cs= " << cs
371 << G4endl;
372 }
373 }
374
375 G4double coeff = 0.0;
376 if(eAverage > 0.) {
377 coeff = cross/eAverage;
378 eAverage /= cross;
379 }
380
381 if(verboseLevel > 1) {
382 G4cout << "Ksi Coefficient for Z= " << Z
383 << " E(keV)= " << lowEdgeEnergy/keV
384 << " Eav(keV)= " << eAverage/keV
385 << " coeff= " << coeff
386 << G4endl;
387 }
388
389 energy->push_back(lowEdgeEnergy);
390 ksi->push_back(coeff);
391 }
392 interp = new G4LogLogInterpolation();
393 G4VEMDataSet* set = new G4EMDataSet(Z,energy,ksi,interp,1.,1.);
394 xsis->AddComponent(set);
395 }
396 if(verboseLevel) xsis->PrintData();
397 shellVacancy->AddXsiTable(xsis);
398 }
399 delete bVector;
400}
401
402
403G4VParticleChange* G4LowEnergyIonisation::PostStepDoIt(const G4Track& track,
404 const G4Step& step)
405{
406 // Delta electron production mechanism on base of the model
407 // J. Stepanek " A program to determine the radiation spectra due
408 // to a single atomic subshell ionisation by a particle or due to
409 // deexcitation or decay of radionuclides",
410 // Comp. Phys. Comm. 1206 pp 1-19 (1997)
411
412 aParticleChange.Initialize(track);
413
414 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
415 G4double kineticEnergy = track.GetKineticEnergy();
416
417 // Select atom and shell
418
419 G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
420 G4int shell = crossSectionHandler->SelectRandomShell(Z, kineticEnergy);
421 const G4AtomicShell* atomicShell =
422 (G4AtomicTransitionManager::Instance())->Shell(Z, shell);
423 G4double bindingEnergy = atomicShell->BindingEnergy();
424 G4int shellId = atomicShell->ShellId();
425
426 // Sample delta energy
427
428 G4int index = couple->GetIndex();
429 G4double tCut = cutForDelta[index];
430 G4double tmax = energySpectrum->MaxEnergyOfSecondaries(kineticEnergy);
431 G4double tDelta = energySpectrum->SampleEnergy(Z, tCut, tmax,
432 kineticEnergy, shell);
433
434 if(tDelta == 0.0)
435 return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
436
437 // Transform to shell potential
438 G4double deltaKinE = tDelta + 2.0*bindingEnergy;
439 G4double primaryKinE = kineticEnergy + 2.0*bindingEnergy;
440
441 // sampling of scattering angle neglecting atomic motion
442 G4double deltaMom = std::sqrt(deltaKinE*(deltaKinE + 2.0*electron_mass_c2));
443 G4double primaryMom = std::sqrt(primaryKinE*(primaryKinE + 2.0*electron_mass_c2));
444
445 G4double cost = deltaKinE * (primaryKinE + 2.0*electron_mass_c2)
446 / (deltaMom * primaryMom);
447
448 if (cost > 1.) cost = 1.;
449 G4double sint = std::sqrt(1. - cost*cost);
450 G4double phi = twopi * G4UniformRand();
451 G4double dirx = sint * std::cos(phi);
452 G4double diry = sint * std::sin(phi);
453 G4double dirz = cost;
454
455 // Rotate to incident electron direction
456 G4ThreeVector primaryDirection = track.GetMomentumDirection();
457 G4ThreeVector deltaDir(dirx,diry,dirz);
458 deltaDir.rotateUz(primaryDirection);
459 dirx = deltaDir.x();
460 diry = deltaDir.y();
461 dirz = deltaDir.z();
462
463
464 // Take into account atomic motion del is relative momentum of the motion
465 // kinetic energy of the motion == bindingEnergy in V.Ivanchenko model
466
467 cost = 2.0*G4UniformRand() - 1.0;
468 sint = std::sqrt(1. - cost*cost);
469 phi = twopi * G4UniformRand();
470 G4double del = std::sqrt(bindingEnergy *(bindingEnergy + 2.0*electron_mass_c2))
471 / deltaMom;
472 dirx += del* sint * std::cos(phi);
473 diry += del* sint * std::sin(phi);
474 dirz += del* cost;
475
476 // Find out new primary electron direction
477 G4double finalPx = primaryMom*primaryDirection.x() - deltaMom*dirx;
478 G4double finalPy = primaryMom*primaryDirection.y() - deltaMom*diry;
479 G4double finalPz = primaryMom*primaryDirection.z() - deltaMom*dirz;
480
481 // create G4DynamicParticle object for delta ray
482 G4DynamicParticle* theDeltaRay = new G4DynamicParticle();
483 theDeltaRay->SetKineticEnergy(tDelta);
484 G4double norm = 1.0/std::sqrt(dirx*dirx + diry*diry + dirz*dirz);
485 dirx *= norm;
486 diry *= norm;
487 dirz *= norm;
488 theDeltaRay->SetMomentumDirection(dirx, diry, dirz);
489 theDeltaRay->SetDefinition(G4Electron::Electron());
490
491 G4double theEnergyDeposit = bindingEnergy;
492
493 // fill ParticleChange
494 // changed energy and momentum of the actual particle
495
496 G4double finalKinEnergy = kineticEnergy - tDelta - theEnergyDeposit;
497 if(finalKinEnergy < 0.0) {
498 theEnergyDeposit += finalKinEnergy;
499 finalKinEnergy = 0.0;
500 aParticleChange.ProposeTrackStatus(fStopAndKill);
501
502 } else {
503
504 G4double norm = 1.0/std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
505 finalPx *= norm;
506 finalPy *= norm;
507 finalPz *= norm;
508 aParticleChange.ProposeMomentumDirection(finalPx, finalPy, finalPz);
509 }
510
511 aParticleChange.ProposeEnergy(finalKinEnergy);
512
513 // Generation of Fluorescence and Auger
514 size_t nSecondaries = 0;
515 size_t totalNumber = 1;
516 std::vector<G4DynamicParticle*>* secondaryVector = 0;
517 G4DynamicParticle* aSecondary = 0;
518 G4ParticleDefinition* type = 0;
519
520 // Fluorescence data start from element 6
521
522 if (Fluorescence() && Z > 5 && (bindingEnergy >= cutForPhotons
523 || bindingEnergy >= cutForElectrons)) {
524
525 secondaryVector = deexcitationManager.GenerateParticles(Z, shellId);
526
527 if (secondaryVector != 0) {
528
529 nSecondaries = secondaryVector->size();
530 for (size_t i = 0; i<nSecondaries; i++) {
531
532 aSecondary = (*secondaryVector)[i];
533 if (aSecondary) {
534
535 G4double e = aSecondary->GetKineticEnergy();
536 type = aSecondary->GetDefinition();
537 if (e < theEnergyDeposit &&
538 ((type == G4Gamma::Gamma() && e > cutForPhotons ) ||
539 (type == G4Electron::Electron() && e > cutForElectrons ))) {
540
541 theEnergyDeposit -= e;
542 totalNumber++;
543
544 } else {
545
546 delete aSecondary;
547 (*secondaryVector)[i] = 0;
548 }
549 }
550 }
551 }
552 }
553
554 // Save delta-electrons
555
556 aParticleChange.SetNumberOfSecondaries(totalNumber);
557 aParticleChange.AddSecondary(theDeltaRay);
558
559 // Save Fluorescence and Auger
560
561 if (secondaryVector) {
562
563 for (size_t l = 0; l < nSecondaries; l++) {
564
565 aSecondary = (*secondaryVector)[l];
566
567 if(aSecondary) {
568
569 aParticleChange.AddSecondary(aSecondary);
570 }
571 }
572 delete secondaryVector;
573 }
574
575 if(theEnergyDeposit < 0.) {
576 G4cout << "G4LowEnergyIonisation: Negative energy deposit: "
577 << theEnergyDeposit/eV << " eV" << G4endl;
578 theEnergyDeposit = 0.0;
579 }
580 aParticleChange.ProposeLocalEnergyDeposit(theEnergyDeposit);
581
582 return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
583}
584
585
586void G4LowEnergyIonisation::PrintInfoDefinition()
587{
588 G4String comments = "Total cross sections from EEDL database.";
589 comments += "\n Gamma energy sampled from a parametrised formula.";
590 comments += "\n Implementation of the continuous dE/dx part.";
591 comments += "\n At present it can be used for electrons ";
592 comments += "in the energy range [250eV,100GeV].";
593 comments += "\n The process must work with G4LowEnergyBremsstrahlung.";
594
595 G4cout << G4endl << GetProcessName() << ": " << comments << G4endl;
596}
597
598G4bool G4LowEnergyIonisation::IsApplicable(const G4ParticleDefinition& particle)
599{
600 return ( (&particle == G4Electron::Electron()) );
601}
602
603std::vector<G4DynamicParticle*>*
604G4LowEnergyIonisation::DeexciteAtom(const G4MaterialCutsCouple* couple,
605 G4double incidentEnergy,
606 G4double eLoss)
607{
608 // create vector of secondary particles
609 const G4Material* material = couple->GetMaterial();
610
611 std::vector<G4DynamicParticle*>* partVector =
612 new std::vector<G4DynamicParticle*>;
613
614 if(eLoss > cutForPhotons && eLoss > cutForElectrons) {
615
616 const G4AtomicTransitionManager* transitionManager =
617 G4AtomicTransitionManager::Instance();
618
619 size_t nElements = material->GetNumberOfElements();
620 const G4ElementVector* theElementVector = material->GetElementVector();
621
622 std::vector<G4DynamicParticle*>* secVector = 0;
623 G4DynamicParticle* aSecondary = 0;
624 G4ParticleDefinition* type = 0;
625 G4double e;
626 G4ThreeVector position;
627 G4int shell, shellId;
628
629 // sample secondaries
630
631 G4double eTot = 0.0;
632 std::vector<G4int> n =
633 shellVacancy->GenerateNumberOfIonisations(couple,
634 incidentEnergy,eLoss);
635 for (size_t i=0; i<nElements; i++) {
636
637 G4int Z = (G4int)((*theElementVector)[i]->GetZ());
638 size_t nVacancies = n[i];
639
640 G4double maxE = transitionManager->Shell(Z, 0)->BindingEnergy();
641
642 if (nVacancies && Z > 5 && (maxE>cutForPhotons || maxE>cutForElectrons)) {
643
644 for (size_t j=0; j<nVacancies; j++) {
645
646 shell = crossSectionHandler->SelectRandomShell(Z, incidentEnergy);
647 shellId = transitionManager->Shell(Z, shell)->ShellId();
648 G4double maxEShell =
649 transitionManager->Shell(Z, shell)->BindingEnergy();
650
651 if (maxEShell>cutForPhotons || maxEShell>cutForElectrons ) {
652
653 secVector = deexcitationManager.GenerateParticles(Z, shellId);
654
655 if (secVector != 0) {
656
657 for (size_t l = 0; l<secVector->size(); l++) {
658
659 aSecondary = (*secVector)[l];
660 if (aSecondary != 0) {
661
662 e = aSecondary->GetKineticEnergy();
663 type = aSecondary->GetDefinition();
664 if ( eTot + e <= eLoss &&
665 (type == G4Gamma::Gamma() && e>cutForPhotons ) ||
666 (type == G4Electron::Electron() && e>cutForElectrons)) {
667
668 eTot += e;
669 partVector->push_back(aSecondary);
670
671 } else {
672
673 delete aSecondary;
674
675 }
676 }
677 }
678 delete secVector;
679 }
680 }
681 }
682 }
683 }
684 }
685 return partVector;
686}
687
688G4double G4LowEnergyIonisation::GetMeanFreePath(const G4Track& track,
689 G4double , // previousStepSize
690 G4ForceCondition* cond)
691{
692 *cond = NotForced;
693 G4int index = (track.GetMaterialCutsCouple())->GetIndex();
694 const G4VEMDataSet* data = theMeanFreePath->GetComponent(index);
695 G4double meanFreePath = data->FindValue(track.GetKineticEnergy());
696 return meanFreePath;
697}
698
699void G4LowEnergyIonisation::SetCutForLowEnSecPhotons(G4double cut)
700{
701 cutForPhotons = cut;
702 deexcitationManager.SetCutForSecondaryPhotons(cut);
703}
704
705void G4LowEnergyIonisation::SetCutForLowEnSecElectrons(G4double cut)
706{
707 cutForElectrons = cut;
708 deexcitationManager.SetCutForAugerElectrons(cut);
709}
710
711void G4LowEnergyIonisation::ActivateAuger(G4bool val)
712{
713 deexcitationManager.ActivateAugerElectronProduction(val);
714}
715
Note: See TracBrowser for help on using the repository browser.