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

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

import all except CVS

File size: 9.7 KB
RevLine 
[819]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: G4LowEnergyCompton.cc,v 1.41 2006/06/29 19:40:15 gunter Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
28//
29// Author: A. Forti
30// Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// History:
33// --------
34// Added Livermore data table construction methods A. Forti
35// Modified BuildMeanFreePath to read new data tables A. Forti
36// Modified PostStepDoIt to insert sampling with EPDL97 data A. Forti
37// Added SelectRandomAtom A. Forti
38// Added map of the elements A. Forti
39// 24.04.2001 V.Ivanchenko - Remove RogueWave
40// 06.08.2001 MGP - Revised according to a design iteration
41// 22.01.2003 V.Ivanchenko - Cut per region
42// 10.03.2003 V.Ivanchenko - Remove CutPerMaterial warning
43// 24.04.2003 V.Ivanchenko - Cut per region mfpt
44//
45// -------------------------------------------------------------------
46
47#include "G4LowEnergyCompton.hh"
48#include "Randomize.hh"
49#include "G4ParticleDefinition.hh"
50#include "G4Track.hh"
51#include "G4Step.hh"
52#include "G4ForceCondition.hh"
53#include "G4Gamma.hh"
54#include "G4Electron.hh"
55#include "G4DynamicParticle.hh"
56#include "G4VParticleChange.hh"
57#include "G4ThreeVector.hh"
58#include "G4EnergyLossTables.hh"
59#include "G4VCrossSectionHandler.hh"
60#include "G4CrossSectionHandler.hh"
61#include "G4VEMDataSet.hh"
62#include "G4CompositeEMDataSet.hh"
63#include "G4VDataSetAlgorithm.hh"
64#include "G4LogLogInterpolation.hh"
65#include "G4VRangeTest.hh"
66#include "G4RangeTest.hh"
67#include "G4MaterialCutsCouple.hh"
68
69
70G4LowEnergyCompton::G4LowEnergyCompton(const G4String& processName)
71 : G4VDiscreteProcess(processName),
72 lowEnergyLimit(250*eV),
73 highEnergyLimit(100*GeV),
74 intrinsicLowEnergyLimit(10*eV),
75 intrinsicHighEnergyLimit(100*GeV)
76{
77 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
78 highEnergyLimit > intrinsicHighEnergyLimit)
79 {
80 G4Exception("G4LowEnergyCompton::G4LowEnergyCompton - energy outside intrinsic process validity range");
81 }
82
83 crossSectionHandler = new G4CrossSectionHandler;
84
85 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
86 G4String scatterFile = "comp/ce-sf-";
87 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
88 scatterFunctionData->LoadData(scatterFile);
89
90 meanFreePathTable = 0;
91
92 rangeTest = new G4RangeTest;
93
94 if (verboseLevel > 0)
95 {
96 G4cout << GetProcessName() << " is created " << G4endl
97 << "Energy range: "
98 << lowEnergyLimit / keV << " keV - "
99 << highEnergyLimit / GeV << " GeV"
100 << G4endl;
101 }
102}
103
104G4LowEnergyCompton::~G4LowEnergyCompton()
105{
106 delete meanFreePathTable;
107 delete crossSectionHandler;
108 delete scatterFunctionData;
109 delete rangeTest;
110}
111
112void G4LowEnergyCompton::BuildPhysicsTable(const G4ParticleDefinition& )
113{
114
115 crossSectionHandler->Clear();
116 G4String crossSectionFile = "comp/ce-cs-";
117 crossSectionHandler->LoadData(crossSectionFile);
118
119 delete meanFreePathTable;
120 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
121}
122
123G4VParticleChange* G4LowEnergyCompton::PostStepDoIt(const G4Track& aTrack,
124 const G4Step& aStep)
125{
126 // The scattered gamma energy is sampled according to Klein - Nishina formula.
127 // then accepted or rejected depending on the Scattering Function multiplied
128 // by factor from Klein - Nishina formula.
129 // Expression of the angular distribution as Klein Nishina
130 // angular and energy distribution and Scattering fuctions is taken from
131 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
132 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
133 // data are interpolated while in the article they are fitted.
134 // Reference to the article is from J. Stepanek New Photon, Positron
135 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
136 // TeV (draft).
137 // The random number techniques of Butcher & Messel are used
138 // (Nucl Phys 20(1960),15).
139
140 aParticleChange.Initialize(aTrack);
141
142 // Dynamic particle quantities
143 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
144 G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
145
146 if (photonEnergy0 <= lowEnergyLimit)
147 {
148 aParticleChange.ProposeTrackStatus(fStopAndKill);
149 aParticleChange.ProposeEnergy(0.);
150 aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
151 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
152 }
153
154 G4double e0m = photonEnergy0 / electron_mass_c2 ;
155 G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
156
157 // Select randomly one element in the current material
158 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
159 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
160
161 G4double epsilon0 = 1. / (1. + 2. * e0m);
162 G4double epsilon0Sq = epsilon0 * epsilon0;
163 G4double alpha1 = -std::log(epsilon0);
164 G4double alpha2 = 0.5 * (1. - epsilon0Sq);
165
166 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
167
168 // Sample the energy of the scattered photon
169 G4double epsilon;
170 G4double epsilonSq;
171 G4double oneCosT;
172 G4double sinT2;
173 G4double gReject;
174 do
175 {
176 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
177 {
178 epsilon = std::exp(-alpha1 * G4UniformRand()); // std::pow(epsilon0,G4UniformRand())
179 epsilonSq = epsilon * epsilon;
180 }
181 else
182 {
183 epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
184 epsilon = std::sqrt(epsilonSq);
185 }
186
187 oneCosT = (1. - epsilon) / ( epsilon * e0m);
188 sinT2 = oneCosT * (2. - oneCosT);
189 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm);
190 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
191 gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
192
193 } while(gReject < G4UniformRand()*Z);
194
195 G4double cosTheta = 1. - oneCosT;
196 G4double sinTheta = std::sqrt (sinT2);
197 G4double phi = twopi * G4UniformRand() ;
198 G4double dirx = sinTheta * std::cos(phi);
199 G4double diry = sinTheta * std::sin(phi);
200 G4double dirz = cosTheta ;
201
202 // Update G4VParticleChange for the scattered photon
203
204 G4ThreeVector photonDirection1(dirx,diry,dirz);
205 photonDirection1.rotateUz(photonDirection0);
206 aParticleChange.ProposeMomentumDirection(photonDirection1) ;
207 G4double photonEnergy1 = epsilon * photonEnergy0;
208
209 if (photonEnergy1 > 0.)
210 {
211 aParticleChange.ProposeEnergy(photonEnergy1) ;
212 }
213 else
214 {
215 aParticleChange.ProposeEnergy(0.) ;
216 aParticleChange.ProposeTrackStatus(fStopAndKill);
217 }
218
219 // Kinematics of the scattered electron
220 G4double eKineticEnergy = photonEnergy0 - photonEnergy1;
221
222 // Generate the electron only if with large enough range w.r.t. cuts and safety
223
224 G4double safety = aStep.GetPostStepPoint()->GetSafety();
225
226 if (rangeTest->Escape(G4Electron::Electron(),couple,eKineticEnergy,safety))
227 {
228 G4double eMomentum = std::sqrt(eKineticEnergy*(eKineticEnergy+2.*electron_mass_c2));
229 G4ThreeVector eDirection((photonEnergy0 * photonDirection0 -
230 photonEnergy1 * photonDirection1) * (1./eMomentum));
231 G4DynamicParticle* electron = new G4DynamicParticle (G4Electron::Electron(),
232 eDirection,eKineticEnergy) ;
233 aParticleChange.SetNumberOfSecondaries(1);
234 aParticleChange.AddSecondary(electron);
235 aParticleChange.ProposeLocalEnergyDeposit(0.);
236 }
237 else
238 {
239 aParticleChange.SetNumberOfSecondaries(0);
240 aParticleChange.ProposeLocalEnergyDeposit(eKineticEnergy);
241 }
242
243 return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
244}
245
246G4bool G4LowEnergyCompton::IsApplicable(const G4ParticleDefinition& particle)
247{
248 return ( &particle == G4Gamma::Gamma() );
249}
250
251G4double G4LowEnergyCompton::GetMeanFreePath(const G4Track& track,
252 G4double, // previousStepSize
253 G4ForceCondition*)
254{
255 const G4DynamicParticle* photon = track.GetDynamicParticle();
256 G4double energy = photon->GetKineticEnergy();
257 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
258 size_t materialIndex = couple->GetIndex();
259
260 G4double meanFreePath;
261 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
262 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
263 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
264 return meanFreePath;
265}
Note: See TracBrowser for help on using the repository browser.