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

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

import all except CVS

File size: 7.6 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//
28// $Id: G4LowEnergyRayleigh.cc,v 1.37 2006/06/29 19:40:29 gunter Exp $
29// GEANT4 tag $Name: $
30//
31// Author: A. Forti
32// Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
33//
34// History:
35// --------
36// Added Livermore data table construction methods A. Forti
37// Added BuildMeanFreePath A. Forti
38// Added PostStepDoIt A. Forti
39// Added SelectRandomAtom A. Forti
40// Added map of the elements A.Forti
41// 24.04.01 V.Ivanchenko remove RogueWave
42// 11.08.2001 MGP - Major revision according to a design iteration
43// 06.10.2001 MGP - Added strategy to test range for secondary generation
44// 05.06.2002 F.Longo and G.Depaola - bug fixed in angular distribution
45// 20.10.2002 G. Depaola - Change sampling method of theta
46// 22.01.2003 V.Ivanchenko - Cut per region
47// 24.04.2003 V.Ivanchenko - Cut per region mfpt
48//
49// --------------------------------------------------------------------
50
51#include "G4LowEnergyRayleigh.hh"
52#include "Randomize.hh"
53#include "G4ParticleDefinition.hh"
54#include "G4Track.hh"
55#include "G4Step.hh"
56#include "G4ForceCondition.hh"
57#include "G4Gamma.hh"
58#include "G4Electron.hh"
59#include "G4DynamicParticle.hh"
60#include "G4VParticleChange.hh"
61#include "G4ThreeVector.hh"
62#include "G4VCrossSectionHandler.hh"
63#include "G4CrossSectionHandler.hh"
64#include "G4VEMDataSet.hh"
65#include "G4CompositeEMDataSet.hh"
66#include "G4VDataSetAlgorithm.hh"
67#include "G4LogLogInterpolation.hh"
68
69#include "G4MaterialCutsCouple.hh"
70
71G4LowEnergyRayleigh::G4LowEnergyRayleigh(const G4String& processName)
72 : G4VDiscreteProcess(processName),
73 lowEnergyLimit(250*eV),
74 highEnergyLimit(100*GeV),
75 intrinsicLowEnergyLimit(10*eV),
76 intrinsicHighEnergyLimit(100*GeV)
77{
78 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
79 highEnergyLimit > intrinsicHighEnergyLimit)
80 {
81 G4Exception("G4LowEnergyRayleigh::G4LowEnergyRayleigh - energy limit outside intrinsic process validity range");
82 }
83
84 crossSectionHandler = new G4CrossSectionHandler();
85
86 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
87 G4String formFactorFile = "rayl/re-ff-";
88 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
89 formFactorData->LoadData(formFactorFile);
90
91 meanFreePathTable = 0;
92
93 if (verboseLevel > 0)
94 {
95 G4cout << GetProcessName() << " is created " << G4endl
96 << "Energy range: "
97 << lowEnergyLimit / keV << " keV - "
98 << highEnergyLimit / GeV << " GeV"
99 << G4endl;
100 }
101}
102
103G4LowEnergyRayleigh::~G4LowEnergyRayleigh()
104{
105 delete meanFreePathTable;
106 delete crossSectionHandler;
107 delete formFactorData;
108}
109
110void G4LowEnergyRayleigh::BuildPhysicsTable(const G4ParticleDefinition& )
111{
112
113 crossSectionHandler->Clear();
114 G4String crossSectionFile = "rayl/re-cs-";
115 crossSectionHandler->LoadData(crossSectionFile);
116
117 delete meanFreePathTable;
118 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
119}
120
121G4VParticleChange* G4LowEnergyRayleigh::PostStepDoIt(const G4Track& aTrack,
122 const G4Step& aStep)
123{
124
125 aParticleChange.Initialize(aTrack);
126
127 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
128 G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
129
130 if (photonEnergy0 <= lowEnergyLimit)
131 {
132 aParticleChange.ProposeTrackStatus(fStopAndKill);
133 aParticleChange.ProposeEnergy(0.);
134 aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
135 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
136 }
137
138 // G4double e0m = photonEnergy0 / electron_mass_c2 ;
139 G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
140
141 // Select randomly one element in the current material
142 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
143 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
144
145 // Sample the angle of the scattered photon
146
147 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
148
149 G4double gReject,x,dataFormFactor;
150 G4double randomFormFactor;
151 G4double cosTheta;
152 G4double sinTheta;
153 G4double fcostheta;
154
155 do
156 {
157 do
158 {
159 cosTheta = 2. * G4UniformRand() - 1.;
160 fcostheta = ( 1. + cosTheta*cosTheta)/2.;
161 } while (fcostheta < G4UniformRand());
162
163 G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.);
164 x = sinThetaHalf / (wlPhoton/cm);
165 if (x > 1.e+005)
166 dataFormFactor = formFactorData->FindValue(x,Z-1);
167 else
168 dataFormFactor = formFactorData->FindValue(0.,Z-1);
169 randomFormFactor = G4UniformRand() * Z * Z;
170 sinTheta = std::sqrt(1. - cosTheta*cosTheta);
171 gReject = dataFormFactor * dataFormFactor;
172
173 } while( gReject < randomFormFactor);
174
175 // Scattered photon angles. ( Z - axis along the parent photon)
176 G4double phi = twopi * G4UniformRand() ;
177 G4double dirX = sinTheta*std::cos(phi);
178 G4double dirY = sinTheta*std::sin(phi);
179 G4double dirZ = cosTheta;
180
181 // Update G4VParticleChange for the scattered photon
182 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
183
184 photonDirection1.rotateUz(photonDirection0);
185 aParticleChange.ProposeEnergy(photonEnergy0);
186 aParticleChange.ProposeMomentumDirection(photonDirection1);
187
188 aParticleChange.SetNumberOfSecondaries(0);
189
190 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
191}
192
193G4bool G4LowEnergyRayleigh::IsApplicable(const G4ParticleDefinition& particle)
194{
195 return ( &particle == G4Gamma::Gamma() );
196}
197
198G4double G4LowEnergyRayleigh::GetMeanFreePath(const G4Track& track,
199 G4double, // previousStepSize
200 G4ForceCondition*)
201{
202 const G4DynamicParticle* photon = track.GetDynamicParticle();
203 G4double energy = photon->GetKineticEnergy();
204 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
205 size_t materialIndex = couple->GetIndex();
206
207 G4double meanFreePath;
208 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
209 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
210 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
211 return meanFreePath;
212}
Note: See TracBrowser for help on using the repository browser.