Geant4  10.01
G4LowEnergyIonisation.cc
Go to the documentation of this file.
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$
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 "G4PhysicalConstants.hh"
107 #include "G4SystemOfUnits.hh"
111 #include "G4RDAtomicShell.hh"
112 #include "G4RDVDataSetAlgorithm.hh"
115 #include "G4RDEMDataSet.hh"
116 #include "G4RDVEMDataSet.hh"
117 #include "G4RDCompositeEMDataSet.hh"
118 #include "G4EnergyLossTables.hh"
119 #include "G4RDShellVacancy.hh"
120 #include "G4UnitsTable.hh"
121 #include "G4Electron.hh"
122 #include "G4Gamma.hh"
123 #include "G4ProductionCutsTable.hh"
124 
126  : G4eLowEnergyLoss(nam),
127  crossSectionHandler(0),
128  theMeanFreePath(0),
129  energySpectrum(0),
130  shellVacancy(0)
131 {
132  cutForPhotons = 250.0*eV;
133  cutForElectrons = 250.0*eV;
134  verboseLevel = 0;
135 }
136 
137 
139 {
140  delete crossSectionHandler;
141  delete energySpectrum;
142  delete theMeanFreePath;
143  delete shellVacancy;
144 }
145 
146 
148 {
149  if(verboseLevel > 0) {
150  G4cout << "G4LowEnergyIonisation::BuildPhysicsTable start"
151  << G4endl;
152  }
153 
154  cutForDelta.clear();
155 
156  // Create and fill IonisationParameters once
157  if( energySpectrum != 0 ) delete energySpectrum;
159 
160  if(verboseLevel > 0) {
161  G4cout << "G4RDVEnergySpectrum is initialized"
162  << G4endl;
163  }
164 
165  // Create and fill G4RDCrossSectionHandler once
166 
167  if ( crossSectionHandler != 0 ) delete crossSectionHandler;
168  G4RDVDataSetAlgorithm* interpolation = new G4RDSemiLogInterpolation();
169  G4double lowKineticEnergy = GetLowerBoundEloss();
170  G4double highKineticEnergy = GetUpperBoundEloss();
171  G4int totBin = GetNbinEloss();
173  interpolation,
174  lowKineticEnergy,
175  highKineticEnergy,
176  totBin);
177  crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
178 
179  if (verboseLevel > 0) {
181  << " is created; Cross section data: "
182  << G4endl;
184  G4cout << "Parameters: "
185  << G4endl;
187  }
188 
189  // Build loss table for IonisationIV
190 
191  BuildLossTable(aParticleType);
192 
193  if(verboseLevel > 0) {
194  G4cout << "The loss table is built"
195  << G4endl;
196  }
197 
198  if (&aParticleType==G4Electron::Electron()) {
199 
200  RecorderOfElectronProcess[CounterOfElectronProcess] = (*this).theLossTable;
203 
204  } else {
205 
206  RecorderOfPositronProcess[CounterOfPositronProcess] = (*this).theLossTable;
208  }
209 
210  // Build mean free path data using cut values
211 
212  if( theMeanFreePath ) delete theMeanFreePath;
214  BuildMeanFreePathForMaterials(&cutForDelta);
215 
216  if(verboseLevel > 0) {
217  G4cout << "The MeanFreePath table is built"
218  << G4endl;
220  }
221 
222  // Build common DEDX table for all ionisation processes
223 
224  BuildDEDXTable(aParticleType);
225 
226  if (verboseLevel > 0) {
227  G4cout << "G4LowEnergyIonisation::BuildPhysicsTable end"
228  << G4endl;
229  }
230 }
231 
232 
234 {
235  // Build table for energy loss due to soft brems
236  // the tables are built for *MATERIALS* binning is taken from LowEnergyLoss
237 
238  G4double lowKineticEnergy = GetLowerBoundEloss();
239  G4double highKineticEnergy = GetUpperBoundEloss();
240  size_t totBin = GetNbinEloss();
241 
242  // create table
243 
244  if (theLossTable) {
246  delete theLossTable;
247  }
248  const G4ProductionCutsTable* theCoupleTable=
250  size_t numOfCouples = theCoupleTable->GetTableSize();
251  theLossTable = new G4PhysicsTable(numOfCouples);
252 
253  if (shellVacancy != 0) delete shellVacancy;
255  G4DataVector* ksi = 0;
256  G4DataVector* energy = 0;
257  size_t binForFluo = totBin/10;
258 
259  G4PhysicsLogVector* bVector = new G4PhysicsLogVector(lowKineticEnergy,
260  highKineticEnergy,
261  binForFluo);
263 
264  // Clean up the vector of cuts
265 
266  cutForDelta.clear();
267 
268  // Loop for materials
269 
270  for (size_t m=0; m<numOfCouples; m++) {
271 
272  // create physics vector and fill it
273  G4PhysicsLogVector* aVector = new G4PhysicsLogVector(lowKineticEnergy,
274  highKineticEnergy,
275  totBin);
276 
277  // get material parameters needed for the energy loss calculation
278  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
279  const G4Material* material= couple->GetMaterial();
280 
281  // the cut cannot be below lowest limit
282  G4double tCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[m];
283  if(tCut > highKineticEnergy) tCut = highKineticEnergy;
284  cutForDelta.push_back(tCut);
285  const G4ElementVector* theElementVector = material->GetElementVector();
286  size_t NumberOfElements = material->GetNumberOfElements() ;
287  const G4double* theAtomicNumDensityVector =
288  material->GetAtomicNumDensityVector();
289  if(verboseLevel > 0) {
290  G4cout << "Energy loss for material # " << m
291  << " tCut(keV)= " << tCut/keV
292  << G4endl;
293  }
294 
295  // now comes the loop for the kinetic energy values
296  for (size_t i = 0; i<totBin; i++) {
297 
298  G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(i);
299  G4double ionloss = 0.;
300 
301  // loop for elements in the material
302  for (size_t iel=0; iel<NumberOfElements; iel++ ) {
303 
304  G4int Z = (G4int)((*theElementVector)[iel]->GetZ());
305 
306  G4int nShells = transitionManager->NumberOfShells(Z);
307 
308  for (G4int n=0; n<nShells; n++) {
309 
310  G4double e = energySpectrum->AverageEnergy(Z, 0.0, tCut,
311  lowEdgeEnergy, n);
312  G4double cs= crossSectionHandler->FindValue(Z, lowEdgeEnergy, n);
313  ionloss += e * cs * theAtomicNumDensityVector[iel];
314 
315  if(verboseLevel > 1 || (Z == 14 && lowEdgeEnergy>1. && lowEdgeEnergy<0.)) {
316  G4cout << "Z= " << Z
317  << " shell= " << n
318  << " E(keV)= " << lowEdgeEnergy/keV
319  << " Eav(keV)= " << e/keV
320  << " cs= " << cs
321  << " loss= " << ionloss
322  << " rho= " << theAtomicNumDensityVector[iel]
323  << G4endl;
324  }
325  }
326  G4double esp = energySpectrum->Excitation(Z, lowEdgeEnergy);
327  ionloss += esp * theAtomicNumDensityVector[iel];
328 
329  }
330  if(verboseLevel > 1 || (m == 0 && lowEdgeEnergy>=1. && lowEdgeEnergy<=0.)) {
331  G4cout << "Sum: "
332  << " E(keV)= " << lowEdgeEnergy/keV
333  << " loss(MeV/mm)= " << ionloss*mm/MeV
334  << G4endl;
335  }
336  aVector->PutValue(i,ionloss);
337  }
338  theLossTable->insert(aVector);
339 
340  // fill data for fluorescence
341 
343  G4RDVEMDataSet* xsis = new G4RDCompositeEMDataSet(interp, 1., 1.);
344  for (size_t iel=0; iel<NumberOfElements; iel++ ) {
345 
346  G4int Z = (G4int)((*theElementVector)[iel]->GetZ());
347  energy = new G4DataVector();
348  ksi = new G4DataVector();
349 
350  for (size_t j = 0; j<binForFluo; j++) {
351 
352  G4double lowEdgeEnergy = bVector->GetLowEdgeEnergy(j);
353  G4double cross = 0.;
354  G4double eAverage= 0.;
355  G4int nShells = transitionManager->NumberOfShells(Z);
356 
357  for (G4int n=0; n<nShells; n++) {
358 
359  G4double e = energySpectrum->AverageEnergy(Z, 0.0, tCut,
360  lowEdgeEnergy, n);
361  G4double pro = energySpectrum->Probability(Z, 0.0, tCut,
362  lowEdgeEnergy, n);
363  G4double cs= crossSectionHandler->FindValue(Z, lowEdgeEnergy, n);
364  eAverage += e * cs * theAtomicNumDensityVector[iel];
365  cross += cs * pro * theAtomicNumDensityVector[iel];
366  if(verboseLevel > 1) {
367  G4cout << "Z= " << Z
368  << " shell= " << n
369  << " E(keV)= " << lowEdgeEnergy/keV
370  << " Eav(keV)= " << e/keV
371  << " pro= " << pro
372  << " cs= " << cs
373  << G4endl;
374  }
375  }
376 
377  G4double coeff = 0.0;
378  if(eAverage > 0.) {
379  coeff = cross/eAverage;
380  eAverage /= cross;
381  }
382 
383  if(verboseLevel > 1) {
384  G4cout << "Ksi Coefficient for Z= " << Z
385  << " E(keV)= " << lowEdgeEnergy/keV
386  << " Eav(keV)= " << eAverage/keV
387  << " coeff= " << coeff
388  << G4endl;
389  }
390 
391  energy->push_back(lowEdgeEnergy);
392  ksi->push_back(coeff);
393  }
394  interp = new G4RDLogLogInterpolation();
395  G4RDVEMDataSet* set = new G4RDEMDataSet(Z,energy,ksi,interp,1.,1.);
396  xsis->AddComponent(set);
397  }
398  if(verboseLevel) xsis->PrintData();
399  shellVacancy->AddXsiTable(xsis);
400  }
401  delete bVector;
402 }
403 
404 
406  const G4Step& step)
407 {
408  // Delta electron production mechanism on base of the model
409  // J. Stepanek " A program to determine the radiation spectra due
410  // to a single atomic subshell ionisation by a particle or due to
411  // deexcitation or decay of radionuclides",
412  // Comp. Phys. Comm. 1206 pp 1-19 (1997)
413 
415 
416  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
417  G4double kineticEnergy = track.GetKineticEnergy();
418 
419  // Select atom and shell
420 
421  G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
422  G4int shell = crossSectionHandler->SelectRandomShell(Z, kineticEnergy);
423  const G4RDAtomicShell* atomicShell =
424  (G4RDAtomicTransitionManager::Instance())->Shell(Z, shell);
425  G4double bindingEnergy = atomicShell->BindingEnergy();
426  G4int shellId = atomicShell->ShellId();
427 
428  // Sample delta energy
429 
430  G4int index = couple->GetIndex();
431  G4double tCut = cutForDelta[index];
432  G4double tmax = energySpectrum->MaxEnergyOfSecondaries(kineticEnergy);
433  G4double tDelta = energySpectrum->SampleEnergy(Z, tCut, tmax,
434  kineticEnergy, shell);
435 
436  if(tDelta == 0.0)
437  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
438 
439  // Transform to shell potential
440  G4double deltaKinE = tDelta + 2.0*bindingEnergy;
441  G4double primaryKinE = kineticEnergy + 2.0*bindingEnergy;
442 
443  // sampling of scattering angle neglecting atomic motion
444  G4double deltaMom = std::sqrt(deltaKinE*(deltaKinE + 2.0*electron_mass_c2));
445  G4double primaryMom = std::sqrt(primaryKinE*(primaryKinE + 2.0*electron_mass_c2));
446 
447  G4double cost = deltaKinE * (primaryKinE + 2.0*electron_mass_c2)
448  / (deltaMom * primaryMom);
449 
450  if (cost > 1.) cost = 1.;
451  G4double sint = std::sqrt(1. - cost*cost);
452  G4double phi = twopi * G4UniformRand();
453  G4double dirx = sint * std::cos(phi);
454  G4double diry = sint * std::sin(phi);
455  G4double dirz = cost;
456 
457  // Rotate to incident electron direction
458  G4ThreeVector primaryDirection = track.GetMomentumDirection();
459  G4ThreeVector deltaDir(dirx,diry,dirz);
460  deltaDir.rotateUz(primaryDirection);
461  dirx = deltaDir.x();
462  diry = deltaDir.y();
463  dirz = deltaDir.z();
464 
465 
466  // Take into account atomic motion del is relative momentum of the motion
467  // kinetic energy of the motion == bindingEnergy in V.Ivanchenko model
468 
469  cost = 2.0*G4UniformRand() - 1.0;
470  sint = std::sqrt(1. - cost*cost);
471  phi = twopi * G4UniformRand();
472  G4double del = std::sqrt(bindingEnergy *(bindingEnergy + 2.0*electron_mass_c2))
473  / deltaMom;
474  dirx += del* sint * std::cos(phi);
475  diry += del* sint * std::sin(phi);
476  dirz += del* cost;
477 
478  // Find out new primary electron direction
479  G4double finalPx = primaryMom*primaryDirection.x() - deltaMom*dirx;
480  G4double finalPy = primaryMom*primaryDirection.y() - deltaMom*diry;
481  G4double finalPz = primaryMom*primaryDirection.z() - deltaMom*dirz;
482 
483  // create G4DynamicParticle object for delta ray
484  G4DynamicParticle* theDeltaRay = new G4DynamicParticle();
485  theDeltaRay->SetKineticEnergy(tDelta);
486  G4double norm = 1.0/std::sqrt(dirx*dirx + diry*diry + dirz*dirz);
487  dirx *= norm;
488  diry *= norm;
489  dirz *= norm;
490  theDeltaRay->SetMomentumDirection(dirx, diry, dirz);
491  theDeltaRay->SetDefinition(G4Electron::Electron());
492 
493  G4double theEnergyDeposit = bindingEnergy;
494 
495  // fill ParticleChange
496  // changed energy and momentum of the actual particle
497 
498  G4double finalKinEnergy = kineticEnergy - tDelta - theEnergyDeposit;
499  if(finalKinEnergy < 0.0) {
500  theEnergyDeposit += finalKinEnergy;
501  finalKinEnergy = 0.0;
503 
504  } else {
505 
506  G4double norm = 1.0/std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
507  finalPx *= norm;
508  finalPy *= norm;
509  finalPz *= norm;
510  aParticleChange.ProposeMomentumDirection(finalPx, finalPy, finalPz);
511  }
512 
513  aParticleChange.ProposeEnergy(finalKinEnergy);
514 
515  // Generation of Fluorescence and Auger
516  size_t nSecondaries = 0;
517  size_t totalNumber = 1;
518  std::vector<G4DynamicParticle*>* secondaryVector = 0;
519  G4DynamicParticle* aSecondary = 0;
520  G4ParticleDefinition* type = 0;
521 
522  // Fluorescence data start from element 6
523 
524  if (Fluorescence() && Z > 5 && (bindingEnergy >= cutForPhotons
525  || bindingEnergy >= cutForElectrons)) {
526 
527  secondaryVector = deexcitationManager.GenerateParticles(Z, shellId);
528 
529  if (secondaryVector != 0) {
530 
531  nSecondaries = secondaryVector->size();
532  for (size_t i = 0; i<nSecondaries; i++) {
533 
534  aSecondary = (*secondaryVector)[i];
535  if (aSecondary) {
536 
537  G4double e = aSecondary->GetKineticEnergy();
538  type = aSecondary->GetDefinition();
539  if (e < theEnergyDeposit &&
540  ((type == G4Gamma::Gamma() && e > cutForPhotons ) ||
541  (type == G4Electron::Electron() && e > cutForElectrons ))) {
542 
543  theEnergyDeposit -= e;
544  totalNumber++;
545 
546  } else {
547 
548  delete aSecondary;
549  (*secondaryVector)[i] = 0;
550  }
551  }
552  }
553  }
554  }
555 
556  // Save delta-electrons
557 
559  aParticleChange.AddSecondary(theDeltaRay);
560 
561  // Save Fluorescence and Auger
562 
563  if (secondaryVector) {
564 
565  for (size_t l = 0; l < nSecondaries; l++) {
566 
567  aSecondary = (*secondaryVector)[l];
568 
569  if(aSecondary) {
570 
571  aParticleChange.AddSecondary(aSecondary);
572  }
573  }
574  delete secondaryVector;
575  }
576 
577  if(theEnergyDeposit < 0.) {
578  G4cout << "G4LowEnergyIonisation: Negative energy deposit: "
579  << theEnergyDeposit/eV << " eV" << G4endl;
580  theEnergyDeposit = 0.0;
581  }
582  aParticleChange.ProposeLocalEnergyDeposit(theEnergyDeposit);
583 
584  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
585 }
586 
587 
589 {
590  G4String comments = "Total cross sections from EEDL database.";
591  comments += "\n Gamma energy sampled from a parametrised formula.";
592  comments += "\n Implementation of the continuous dE/dx part.";
593  comments += "\n At present it can be used for electrons ";
594  comments += "in the energy range [250eV,100GeV].";
595  comments += "\n The process must work with G4LowEnergyBremsstrahlung.";
596 
597  G4cout << G4endl << GetProcessName() << ": " << comments << G4endl;
598 }
599 
601 {
602  return ( (&particle == G4Electron::Electron()) );
603 }
604 
605 std::vector<G4DynamicParticle*>*
607  G4double incidentEnergy,
608  G4double eLoss)
609 {
610  // create vector of secondary particles
611  const G4Material* material = couple->GetMaterial();
612 
613  std::vector<G4DynamicParticle*>* partVector =
614  new std::vector<G4DynamicParticle*>;
615 
616  if(eLoss > cutForPhotons && eLoss > cutForElectrons) {
617 
618  const G4RDAtomicTransitionManager* transitionManager =
620 
621  size_t nElements = material->GetNumberOfElements();
622  const G4ElementVector* theElementVector = material->GetElementVector();
623 
624  std::vector<G4DynamicParticle*>* secVector = 0;
625  G4DynamicParticle* aSecondary = 0;
626  G4ParticleDefinition* type = 0;
627  G4double e;
629  G4int shell, shellId;
630 
631  // sample secondaries
632 
633  G4double eTot = 0.0;
634  std::vector<G4int> n =
636  incidentEnergy,eLoss);
637  for (size_t i=0; i<nElements; i++) {
638 
639  G4int Z = (G4int)((*theElementVector)[i]->GetZ());
640  size_t nVacancies = n[i];
641 
642  G4double maxE = transitionManager->Shell(Z, 0)->BindingEnergy();
643 
644  if (nVacancies && Z > 5 && (maxE>cutForPhotons || maxE>cutForElectrons)) {
645 
646  for (size_t j=0; j<nVacancies; j++) {
647 
648  shell = crossSectionHandler->SelectRandomShell(Z, incidentEnergy);
649  shellId = transitionManager->Shell(Z, shell)->ShellId();
650  G4double maxEShell =
651  transitionManager->Shell(Z, shell)->BindingEnergy();
652 
653  if (maxEShell>cutForPhotons || maxEShell>cutForElectrons ) {
654 
655  secVector = deexcitationManager.GenerateParticles(Z, shellId);
656 
657  if (secVector != 0) {
658 
659  for (size_t l = 0; l<secVector->size(); l++) {
660 
661  aSecondary = (*secVector)[l];
662  if (aSecondary != 0) {
663 
664  e = aSecondary->GetKineticEnergy();
665  type = aSecondary->GetDefinition();
666  if ( eTot + e <= eLoss &&
667  ((type == G4Gamma::Gamma() && e>cutForPhotons ) ||
668  (type == G4Electron::Electron() && e>cutForElectrons))) {
669 
670  eTot += e;
671  partVector->push_back(aSecondary);
672 
673  } else {
674 
675  delete aSecondary;
676 
677  }
678  }
679  }
680  delete secVector;
681  }
682  }
683  }
684  }
685  }
686  }
687  return partVector;
688 }
689 
691  G4double , // previousStepSize
692  G4ForceCondition* cond)
693 {
694  *cond = NotForced;
695  G4int index = (track.GetMaterialCutsCouple())->GetIndex();
696  const G4RDVEMDataSet* data = theMeanFreePath->GetComponent(index);
697  G4double meanFreePath = data->FindValue(track.GetKineticEnergy());
698  return meanFreePath;
699 }
700 
702 {
703  cutForPhotons = cut;
705 }
706 
708 {
709  cutForElectrons = cut;
711 }
712 
714 {
716 }
717 
G4RDShellVacancy * shellVacancy
static const double MeV
Definition: G4SIunits.hh:193
virtual G4double Excitation(G4int Z, G4double kineticEnergy) const =0
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
void AddXsiTable(G4RDVEMDataSet *set)
void SetCutForSecondaryPhotons(G4double cut)
G4int verboseLevel
Definition: G4VProcess.hh:368
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
std::vector< G4Element * > G4ElementVector
void insert(G4PhysicsVector *)
G4RDVEnergySpectrum * energySpectrum
void ActivateAugerElectronProduction(G4bool val)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
static G4int CounterOfElectronProcess
G4bool Fluorescence() const
virtual std::vector< G4DynamicParticle * > * DeexciteAtom(const G4MaterialCutsCouple *couple, G4double incidentEnergy, G4double eLoss)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void LoadShellData(const G4String &dataFile)
virtual const G4RDVEMDataSet * GetComponent(G4int componentId) const =0
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4ParticleDefinition * GetDefinition() const
static G4double GetLowerBoundEloss()
G4double FindValue(G4int Z, G4double e) const
G4double GetLowEdgeEnergy(size_t binNumber) const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
void BuildPhysicsTable(const G4ParticleDefinition &ParticleType)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static G4RDAtomicTransitionManager * Instance()
std::vector< G4int > GenerateNumberOfIonisations(const G4MaterialCutsCouple *couple, G4double incidentEnergy, G4double eLoss) const
G4double GetKineticEnergy() const
#define position
Definition: xmlparse.cc:605
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
G4bool IsApplicable(const G4ParticleDefinition &)
void SetCutForLowEnSecPhotons(G4double cut)
bool G4bool
Definition: G4Types.hh:79
G4double BindingEnergy() const
G4int SelectRandomShell(G4int Z, G4double e) const
void PutValue(size_t index, G4double theValue)
void BuildDEDXTable(const G4ParticleDefinition &aParticleType)
G4LowEnergyIonisation(const G4String &processName="LowEnergyIoni")
Definition: G4Step.hh:76
const G4int n
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual void PrintData(void) const =0
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual void PrintData() const =0
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
void SetKineticEnergy(G4double aEnergy)
G4RDAtomicDeexcitation deexcitationManager
void SetCutForLowEnSecElectrons(G4double cut)
static G4PhysicsTable ** RecorderOfPositronProcess
virtual void Initialize(const G4Track &)
static G4ProductionCutsTable * GetProductionCutsTable()
static G4double GetUpperBoundEloss()
virtual G4double AverageEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
static const double eV
Definition: G4SIunits.hh:194
static G4PhysicsTable ** RecorderOfElectronProcess
const G4ThreeVector & GetMomentumDirection() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
G4double energy(const ThreeVector &p, const G4double m)
void SetNumberOfSecondaries(G4int totSecondaries)
static G4int GetNbinEloss()
void SetCutForAugerElectrons(G4double cut)
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4RDAtomicShell * Shell(G4int Z, size_t shellIndex) const
void AddSecondary(G4Track *aSecondary)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:110
G4RDVEMDataSet * theMeanFreePath
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static const double keV
Definition: G4SIunits.hh:195
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4PhysicsTable * theLossTable
static G4int CounterOfPositronProcess
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ForceCondition
virtual void AddComponent(G4RDVEMDataSet *dataSet)=0
G4double bindingEnergy(G4int A, G4int Z)
void BuildLossTable(const G4ParticleDefinition &ParticleType)
static const double mm
Definition: G4SIunits.hh:102
G4RDVCrossSectionHandler * crossSectionHandler
virtual G4double SampleEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
void clearAndDestroy()
const G4Material * GetMaterial() const
G4int ShellId() const
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
virtual G4double MaxEnergyOfSecondaries(G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=0) const =0