source: trunk/source/processes/electromagnetic/lowenergy/src/G4LowEnergyPhotoElectric.cc@ 966

Last change on this file since 966 was 961, checked in by garnier, 17 years ago

update processes

File size: 14.3 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//
27// $Id: G4LowEnergyPhotoElectric.cc,v 1.56 2006/06/29 19:40:23 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// Author: A. Forti
31// Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
32//
33// History:
34// --------
35// October 1998 - low energy modifications by Alessandra Forti
36// Added Livermore data table construction methods A. Forti
37// Modified BuildMeanFreePath to read new data tables A. Forti
38// Added EnergySampling method A. Forti
39// Modified PostStepDoIt to insert sampling with EPDL97 data A. Forti
40// Added SelectRandomAtom A. Forti
41// Added map of the elements A. Forti
42// 10.04.2000 VL
43// - Correcting Fluorescence transition probabilities in order to take into account
44// non-radiative transitions. No Auger electron simulated yet: energy is locally deposited.
45// 17.02.2000 Veronique Lefebure
46// - bugs corrected in fluorescence simulation:
47// . when final use of binding energy: no photon was ever created
48// . no Fluorescence was simulated when the photo-electron energy
49// was below production threshold.
50//
51// 07-09-99, if no e- emitted: edep=photon energy, mma
52// 24.04.01 V.Ivanchenko remove RogueWave
53// 12.08.2001 MGP Revised according to a design iteration
54// 16.09.2001 E. Guardincerri Added fluorescence generation
55// 06.10.2001 MGP Added protection to avoid negative electron energies
56// when binding energy of selected shell > photon energy
57// 18.04.2001 V.Ivanchenko Fix problem with low energy gammas from fluorescence
58// MeanFreePath is calculated by crosSectionHandler directly
59// 31.05.2002 V.Ivanchenko Add path of Fluo + Auger cuts to AtomicDeexcitation
60// 14.06.2002 V.Ivanchenko By default do not cheak range of e-
61// 21.01.2003 V.Ivanchenko Cut per region
62// 10.05.2004 P.Rodrigues Changes to accommodate new angular generators
63// 20.01.2006 A.Trindade Changes to accommodate polarized angular generator
64//
65// --------------------------------------------------------------
66
67#include "G4LowEnergyPhotoElectric.hh"
68
69#include "G4VPhotoElectricAngularDistribution.hh"
70#include "G4PhotoElectricAngularGeneratorSimple.hh"
71#include "G4PhotoElectricAngularGeneratorSauterGavrila.hh"
72#include "G4PhotoElectricAngularGeneratorPolarized.hh"
73
74#include "G4ParticleDefinition.hh"
75#include "G4Track.hh"
76#include "G4Step.hh"
77#include "G4ForceCondition.hh"
78#include "G4Gamma.hh"
79#include "G4Electron.hh"
80#include "G4DynamicParticle.hh"
81#include "G4VParticleChange.hh"
82#include "G4ThreeVector.hh"
83#include "G4VCrossSectionHandler.hh"
84#include "G4CrossSectionHandler.hh"
85#include "G4VEMDataSet.hh"
86#include "G4CompositeEMDataSet.hh"
87#include "G4VDataSetAlgorithm.hh"
88#include "G4LogLogInterpolation.hh"
89#include "G4VRangeTest.hh"
90#include "G4RangeNoTest.hh"
91#include "G4AtomicTransitionManager.hh"
92#include "G4AtomicShell.hh"
93#include "G4ProductionCutsTable.hh"
94
95G4LowEnergyPhotoElectric::G4LowEnergyPhotoElectric(const G4String& processName)
96 : G4VDiscreteProcess(processName), lowEnergyLimit(250*eV), highEnergyLimit(100*GeV),
97 intrinsicLowEnergyLimit(10*eV),
98 intrinsicHighEnergyLimit(100*GeV),
99 cutForLowEnergySecondaryPhotons(250.*eV),
100 cutForLowEnergySecondaryElectrons(250.*eV)
101{
102 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
103 highEnergyLimit > intrinsicHighEnergyLimit)
104 {
105 G4Exception("G4LowEnergyPhotoElectric::G4LowEnergyPhotoElectric - energy limit outside intrinsic process validity range");
106 }
107
108 crossSectionHandler = new G4CrossSectionHandler();
109 shellCrossSectionHandler = new G4CrossSectionHandler();
110 meanFreePathTable = 0;
111 rangeTest = new G4RangeNoTest;
112 generatorName = "geant4.6.2";
113 ElectronAngularGenerator = new G4PhotoElectricAngularGeneratorSimple("GEANTSimpleGenerator"); // default generator
114
115
116 if (verboseLevel > 0)
117 {
118 G4cout << GetProcessName() << " is created " << G4endl
119 << "Energy range: "
120 << lowEnergyLimit / keV << " keV - "
121 << highEnergyLimit / GeV << " GeV"
122 << G4endl;
123 }
124}
125
126G4LowEnergyPhotoElectric::~G4LowEnergyPhotoElectric()
127{
128 delete crossSectionHandler;
129 delete shellCrossSectionHandler;
130 delete meanFreePathTable;
131 delete rangeTest;
132 delete ElectronAngularGenerator;
133}
134
135void G4LowEnergyPhotoElectric::BuildPhysicsTable(const G4ParticleDefinition& )
136{
137
138 crossSectionHandler->Clear();
139 G4String crossSectionFile = "phot/pe-cs-";
140 crossSectionHandler->LoadData(crossSectionFile);
141
142 shellCrossSectionHandler->Clear();
143 G4String shellCrossSectionFile = "phot/pe-ss-cs-";
144 shellCrossSectionHandler->LoadShellData(shellCrossSectionFile);
145
146 delete meanFreePathTable;
147 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
148}
149
150G4VParticleChange* G4LowEnergyPhotoElectric::PostStepDoIt(const G4Track& aTrack,
151 const G4Step& aStep)
152{
153 // Fluorescence generated according to:
154 // J. Stepanek ,"A program to determine the radiation spectra due to a single atomic
155 // subshell ionisation by a particle or due to deexcitation or decay of radionuclides",
156 // Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
157
158 aParticleChange.Initialize(aTrack);
159
160 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
161 G4double photonEnergy = incidentPhoton->GetKineticEnergy();
162 if (photonEnergy <= lowEnergyLimit)
163 {
164 aParticleChange.ProposeTrackStatus(fStopAndKill);
165 aParticleChange.ProposeEnergy(0.);
166 aParticleChange.ProposeLocalEnergyDeposit(photonEnergy);
167 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
168 }
169
170 G4ThreeVector photonDirection = incidentPhoton->GetMomentumDirection(); // Returns the normalized direction of the momentum
171
172 // Select randomly one element in the current material
173 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
174 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy);
175
176 // Select the ionised shell in the current atom according to shell cross sections
177 size_t shellIndex = shellCrossSectionHandler->SelectRandomShell(Z,photonEnergy);
178
179 // Retrieve the corresponding identifier and binding energy of the selected shell
180 const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
181 const G4AtomicShell* shell = transitionManager->Shell(Z,shellIndex);
182 G4double bindingEnergy = shell->BindingEnergy();
183 G4int shellId = shell->ShellId();
184
185 // Create lists of pointers to DynamicParticles (photons and electrons)
186 // (Is the electron vector necessary? To be checked)
187 std::vector<G4DynamicParticle*>* photonVector = 0;
188 std::vector<G4DynamicParticle*> electronVector;
189
190 G4double energyDeposit = 0.0;
191
192 // Primary outcoming electron
193 G4double eKineticEnergy = photonEnergy - bindingEnergy;
194
195 // There may be cases where the binding energy of the selected shell is > photon energy
196 // In such cases do not generate secondaries
197 if (eKineticEnergy > 0.)
198 {
199 // Generate the electron only if with large enough range w.r.t. cuts and safety
200 G4double safety = aStep.GetPostStepPoint()->GetSafety();
201
202 if (rangeTest->Escape(G4Electron::Electron(),couple,eKineticEnergy,safety))
203 {
204
205 // Calculate direction of the photoelectron
206 G4ThreeVector gammaPolarization = incidentPhoton->GetPolarization();
207 G4ThreeVector electronDirection = ElectronAngularGenerator->GetPhotoElectronDirection(photonDirection,eKineticEnergy,gammaPolarization,shellId);
208
209 // The electron is created ...
210 G4DynamicParticle* electron = new G4DynamicParticle (G4Electron::Electron(),
211 electronDirection,
212 eKineticEnergy);
213 electronVector.push_back(electron);
214 }
215 else
216 {
217 energyDeposit += eKineticEnergy;
218 }
219 }
220 else
221 {
222 bindingEnergy = photonEnergy;
223 }
224
225 G4int nElectrons = electronVector.size();
226 size_t nTotPhotons = 0;
227 G4int nPhotons=0;
228 const G4ProductionCutsTable* theCoupleTable=
229 G4ProductionCutsTable::GetProductionCutsTable();
230
231 size_t index = couple->GetIndex();
232 G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[index];
233 cutg = std::min(cutForLowEnergySecondaryPhotons,cutg);
234
235 G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[index];
236 cute = std::min(cutForLowEnergySecondaryPhotons,cute);
237
238 G4DynamicParticle* aPhoton;
239
240 // Generation of fluorescence
241 // Data in EADL are available only for Z > 5
242 // Protection to avoid generating photons in the unphysical case of
243 // shell binding energy > photon energy
244 if (Z > 5 && (bindingEnergy > cutg || bindingEnergy > cute))
245 {
246 photonVector = deexcitationManager.GenerateParticles(Z,shellId);
247 nTotPhotons = photonVector->size();
248 for (size_t k=0; k<nTotPhotons; k++)
249 {
250 aPhoton = (*photonVector)[k];
251 if (aPhoton)
252 {
253 G4double itsCut = cutg;
254 if(aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cute;
255 G4double itsEnergy = aPhoton->GetKineticEnergy();
256
257 if (itsEnergy > itsCut && itsEnergy <= bindingEnergy)
258 {
259 nPhotons++;
260 // Local energy deposit is given as the sum of the
261 // energies of incident photons minus the energies
262 // of the outcoming fluorescence photons
263 bindingEnergy -= itsEnergy;
264
265 }
266 else
267 {
268 delete aPhoton;
269 (*photonVector)[k] = 0;
270 }
271 }
272 }
273 }
274
275 energyDeposit += bindingEnergy;
276
277 G4int nSecondaries = nElectrons + nPhotons;
278 aParticleChange.SetNumberOfSecondaries(nSecondaries);
279
280 for (G4int l = 0; l<nElectrons; l++ )
281 {
282 aPhoton = electronVector[l];
283 if(aPhoton) {
284 aParticleChange.AddSecondary(aPhoton);
285 }
286 }
287 for ( size_t ll = 0; ll < nTotPhotons; ll++)
288 {
289 aPhoton = (*photonVector)[ll];
290 if(aPhoton) {
291 aParticleChange.AddSecondary(aPhoton);
292 }
293 }
294
295 delete photonVector;
296
297 if (energyDeposit < 0)
298 {
299 G4cout << "WARNING - "
300 << "G4LowEnergyPhotoElectric::PostStepDoIt - Negative energy deposit"
301 << G4endl;
302 energyDeposit = 0;
303 }
304
305 // Kill the incident photon
306 aParticleChange.ProposeMomentumDirection( 0., 0., 0. );
307 aParticleChange.ProposeEnergy( 0. );
308
309 aParticleChange.ProposeLocalEnergyDeposit(energyDeposit);
310 aParticleChange.ProposeTrackStatus( fStopAndKill );
311
312 // Reset NbOfInteractionLengthLeft and return aParticleChange
313 return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
314}
315
316G4bool G4LowEnergyPhotoElectric::IsApplicable(const G4ParticleDefinition& particle)
317{
318 return ( &particle == G4Gamma::Gamma() );
319}
320
321G4double G4LowEnergyPhotoElectric::GetMeanFreePath(const G4Track& track,
322 G4double, // previousStepSize
323 G4ForceCondition*)
324{
325 const G4DynamicParticle* photon = track.GetDynamicParticle();
326 G4double energy = photon->GetKineticEnergy();
327 G4Material* material = track.GetMaterial();
328 // size_t materialIndex = material->GetIndex();
329
330 G4double meanFreePath = DBL_MAX;
331
332 // if (energy > highEnergyLimit)
333 // meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
334 // else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
335 // else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
336
337 G4double cross = shellCrossSectionHandler->ValueForMaterial(material,energy);
338 if(cross > 0.0) meanFreePath = 1.0/cross;
339
340 return meanFreePath;
341}
342
343void G4LowEnergyPhotoElectric::SetCutForLowEnSecPhotons(G4double cut)
344{
345 cutForLowEnergySecondaryPhotons = cut;
346 deexcitationManager.SetCutForSecondaryPhotons(cut);
347}
348
349void G4LowEnergyPhotoElectric::SetCutForLowEnSecElectrons(G4double cut)
350{
351 cutForLowEnergySecondaryElectrons = cut;
352 deexcitationManager.SetCutForAugerElectrons(cut);
353}
354
355void G4LowEnergyPhotoElectric::ActivateAuger(G4bool val)
356{
357 deexcitationManager.ActivateAugerElectronProduction(val);
358}
359
360void G4LowEnergyPhotoElectric::SetAngularGenerator(G4VPhotoElectricAngularDistribution* distribution)
361{
362 ElectronAngularGenerator = distribution;
363 ElectronAngularGenerator->PrintGeneratorInformation();
364}
365
366void G4LowEnergyPhotoElectric::SetAngularGenerator(const G4String& name)
367{
368 if (name == "default")
369 {
370 delete ElectronAngularGenerator;
371 ElectronAngularGenerator = new G4PhotoElectricAngularGeneratorSimple("GEANT4LowEnergySimpleGenerator");
372 generatorName = name;
373 }
374 else if (name == "standard")
375 {
376 delete ElectronAngularGenerator;
377 ElectronAngularGenerator = new G4PhotoElectricAngularGeneratorSauterGavrila("GEANT4SauterGavrilaGenerator");
378 generatorName = name;
379 }
380 else if (name == "polarized")
381 {
382 delete ElectronAngularGenerator;
383 ElectronAngularGenerator = new G4PhotoElectricAngularGeneratorPolarized("GEANT4LowEnergyPolarizedGenerator");
384 generatorName = name;
385 }
386 else
387 {
388 G4Exception("G4LowEnergyPhotoElectric::SetAngularGenerator - generator does not exist");
389 }
390
391 ElectronAngularGenerator->PrintGeneratorInformation();
392}
Note: See TracBrowser for help on using the repository browser.