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: G4LivermoreRayleighModel.cc,v 1.2 2009/01/21 10:58:13 sincerti Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
28 | // |
---|
29 | |
---|
30 | #include "G4LivermoreRayleighModel.hh" |
---|
31 | |
---|
32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
33 | |
---|
34 | using namespace std; |
---|
35 | |
---|
36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
37 | |
---|
38 | G4LivermoreRayleighModel::G4LivermoreRayleighModel(const G4ParticleDefinition*, |
---|
39 | const G4String& nam) |
---|
40 | :G4VEmModel(nam),isInitialised(false),meanFreePathTable(0),formFactorData(0),crossSectionHandler(0) |
---|
41 | { |
---|
42 | lowEnergyLimit = 250 * eV; // SI - Could be 10 eV ? |
---|
43 | highEnergyLimit = 100 * GeV; |
---|
44 | |
---|
45 | SetLowEnergyLimit(lowEnergyLimit); |
---|
46 | SetHighEnergyLimit(highEnergyLimit); |
---|
47 | // |
---|
48 | verboseLevel= 0; |
---|
49 | // Verbosity scale: |
---|
50 | // 0 = nothing |
---|
51 | // 1 = warning for energy non-conservation |
---|
52 | // 2 = details of energy budget |
---|
53 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
54 | // 4 = entering in methods |
---|
55 | |
---|
56 | G4cout << "Livermore Rayleigh is constructed " << G4endl |
---|
57 | << "Energy range: " |
---|
58 | << lowEnergyLimit / keV << " keV - " |
---|
59 | << highEnergyLimit / GeV << " GeV" |
---|
60 | << G4endl; |
---|
61 | } |
---|
62 | |
---|
63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
64 | |
---|
65 | G4LivermoreRayleighModel::~G4LivermoreRayleighModel() |
---|
66 | { |
---|
67 | if (meanFreePathTable) delete meanFreePathTable; |
---|
68 | if (crossSectionHandler) delete crossSectionHandler; |
---|
69 | if (formFactorData) delete formFactorData; |
---|
70 | } |
---|
71 | |
---|
72 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
73 | |
---|
74 | void G4LivermoreRayleighModel::Initialise(const G4ParticleDefinition* particle, |
---|
75 | const G4DataVector& cuts) |
---|
76 | { |
---|
77 | if (verboseLevel > 3) |
---|
78 | G4cout << "Calling G4LivermoreRayleighModel::Initialise()" << G4endl; |
---|
79 | |
---|
80 | if (crossSectionHandler) |
---|
81 | { |
---|
82 | crossSectionHandler->Clear(); |
---|
83 | delete crossSectionHandler; |
---|
84 | } |
---|
85 | |
---|
86 | // Energy limits |
---|
87 | |
---|
88 | if (LowEnergyLimit() < lowEnergyLimit) |
---|
89 | { |
---|
90 | G4cout << "G4LivermoreRayleighModel: low energy limit increased from " << |
---|
91 | LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl; |
---|
92 | SetLowEnergyLimit(lowEnergyLimit); |
---|
93 | } |
---|
94 | |
---|
95 | if (HighEnergyLimit() > highEnergyLimit) |
---|
96 | { |
---|
97 | G4cout << "G4LivermoreRayleighModel: high energy limit decreased from " << |
---|
98 | HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl; |
---|
99 | SetHighEnergyLimit(highEnergyLimit); |
---|
100 | } |
---|
101 | |
---|
102 | // Data are read for all materials |
---|
103 | |
---|
104 | crossSectionHandler = new G4CrossSectionHandler; |
---|
105 | crossSectionHandler->Clear(); |
---|
106 | G4String crossSectionFile = "rayl/re-cs-"; |
---|
107 | crossSectionHandler->LoadData(crossSectionFile); |
---|
108 | |
---|
109 | meanFreePathTable = 0; |
---|
110 | meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials(); |
---|
111 | |
---|
112 | G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation; |
---|
113 | G4String formFactorFile = "rayl/re-ff-"; |
---|
114 | formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.); |
---|
115 | formFactorData->LoadData(formFactorFile); |
---|
116 | |
---|
117 | // |
---|
118 | |
---|
119 | if (verboseLevel > 2) |
---|
120 | G4cout << "Loaded cross section files for Livermore Rayleigh model" << G4endl; |
---|
121 | |
---|
122 | InitialiseElementSelectors(particle,cuts); |
---|
123 | |
---|
124 | G4cout << "Livermore Rayleigh model is initialized " << G4endl |
---|
125 | << "Energy range: " |
---|
126 | << LowEnergyLimit() / keV << " keV - " |
---|
127 | << HighEnergyLimit() / GeV << " GeV" |
---|
128 | << G4endl; |
---|
129 | |
---|
130 | if(isInitialised) return; |
---|
131 | |
---|
132 | if(pParticleChange) |
---|
133 | fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
134 | else |
---|
135 | fParticleChange = new G4ParticleChangeForGamma(); |
---|
136 | |
---|
137 | isInitialised = true; |
---|
138 | |
---|
139 | } |
---|
140 | |
---|
141 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
142 | |
---|
143 | G4double G4LivermoreRayleighModel::ComputeCrossSectionPerAtom( |
---|
144 | const G4ParticleDefinition*, |
---|
145 | G4double GammaEnergy, |
---|
146 | G4double Z, G4double, |
---|
147 | G4double, G4double) |
---|
148 | { |
---|
149 | if (verboseLevel > 3) |
---|
150 | G4cout << "Calling CrossSectionPerAtom() of G4LivermoreRayleighModel" << G4endl; |
---|
151 | |
---|
152 | G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy); |
---|
153 | return cs; |
---|
154 | } |
---|
155 | |
---|
156 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
157 | |
---|
158 | void G4LivermoreRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/, |
---|
159 | const G4MaterialCutsCouple* couple, |
---|
160 | const G4DynamicParticle* aDynamicGamma, |
---|
161 | G4double, |
---|
162 | G4double) |
---|
163 | { |
---|
164 | if (verboseLevel > 3) |
---|
165 | G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel" << G4endl; |
---|
166 | |
---|
167 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
168 | |
---|
169 | if (photonEnergy0 <= lowEnergyLimit) |
---|
170 | { |
---|
171 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
172 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
173 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); |
---|
174 | return ; |
---|
175 | } |
---|
176 | |
---|
177 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
178 | |
---|
179 | // Select randomly one element in the current material |
---|
180 | G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0); |
---|
181 | |
---|
182 | // Sample the angle of the scattered photon |
---|
183 | |
---|
184 | G4double wlPhoton = h_Planck*c_light/photonEnergy0; |
---|
185 | |
---|
186 | G4double gReject,x,dataFormFactor; |
---|
187 | G4double randomFormFactor; |
---|
188 | G4double cosTheta; |
---|
189 | G4double sinTheta; |
---|
190 | G4double fcostheta; |
---|
191 | |
---|
192 | do |
---|
193 | { |
---|
194 | do |
---|
195 | { |
---|
196 | cosTheta = 2. * G4UniformRand() - 1.; |
---|
197 | fcostheta = ( 1. + cosTheta*cosTheta)/2.; |
---|
198 | } while (fcostheta < G4UniformRand()); |
---|
199 | |
---|
200 | G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.); |
---|
201 | x = sinThetaHalf / (wlPhoton/cm); |
---|
202 | if (x > 1.e+005) |
---|
203 | dataFormFactor = formFactorData->FindValue(x,Z-1); |
---|
204 | else |
---|
205 | dataFormFactor = formFactorData->FindValue(0.,Z-1); |
---|
206 | randomFormFactor = G4UniformRand() * Z * Z; |
---|
207 | sinTheta = std::sqrt(1. - cosTheta*cosTheta); |
---|
208 | gReject = dataFormFactor * dataFormFactor; |
---|
209 | |
---|
210 | } while( gReject < randomFormFactor); |
---|
211 | |
---|
212 | // Scattered photon angles. ( Z - axis along the parent photon) |
---|
213 | G4double phi = twopi * G4UniformRand() ; |
---|
214 | G4double dirX = sinTheta*std::cos(phi); |
---|
215 | G4double dirY = sinTheta*std::sin(phi); |
---|
216 | G4double dirZ = cosTheta; |
---|
217 | |
---|
218 | // Update G4VParticleChange for the scattered photon |
---|
219 | G4ThreeVector photonDirection1(dirX, dirY, dirZ); |
---|
220 | photonDirection1.rotateUz(photonDirection0); |
---|
221 | fParticleChange->ProposeMomentumDirection(photonDirection1); |
---|
222 | |
---|
223 | fParticleChange->SetProposedKineticEnergy(photonEnergy0); |
---|
224 | } |
---|
225 | |
---|
226 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
227 | |
---|
228 | G4double G4LivermoreRayleighModel::GetMeanFreePath(const G4Track& track, |
---|
229 | G4double, // previousStepSize |
---|
230 | G4ForceCondition*) |
---|
231 | { |
---|
232 | const G4DynamicParticle* photon = track.GetDynamicParticle(); |
---|
233 | G4double energy = photon->GetKineticEnergy(); |
---|
234 | const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple(); |
---|
235 | size_t materialIndex = couple->GetIndex(); |
---|
236 | |
---|
237 | G4double meanFreePath; |
---|
238 | if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex); |
---|
239 | else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX; |
---|
240 | else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex); |
---|
241 | return meanFreePath; |
---|
242 | } |
---|
243 | |
---|