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: G4LowEnergyPolarizedRayleigh.cc,v 1.10 2009/06/11 15:47:08 mantero Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
28 | // |
---|
29 | // -------------------------------------------------------------- |
---|
30 | // |
---|
31 | // File name: G4LowEnergyPolarizedRayleigh.cc |
---|
32 | // |
---|
33 | // Author: Capra Riccardo |
---|
34 | // |
---|
35 | // Creation date: May 2005 |
---|
36 | // |
---|
37 | // History: |
---|
38 | // ----------- |
---|
39 | // 02 May 2005 R. Capra 1st implementation |
---|
40 | // 20 May 2005 MGP Changed name of a local variable hiding |
---|
41 | // a data member of the base class |
---|
42 | // |
---|
43 | //---------------------------------------------------------------- |
---|
44 | |
---|
45 | |
---|
46 | #include "G4LowEnergyPolarizedRayleigh.hh" |
---|
47 | |
---|
48 | #include "G4LogLogInterpolation.hh" |
---|
49 | #include "G4VCrossSectionHandler.hh" |
---|
50 | #include "G4VEMDataSet.hh" |
---|
51 | #include "globals.hh" //G4Exception |
---|
52 | |
---|
53 | |
---|
54 | G4LowEnergyPolarizedRayleigh::G4LowEnergyPolarizedRayleigh(const G4String& processName) |
---|
55 | : |
---|
56 | G4VLowEnergyDiscretePhotonProcess(processName, "rayl/re-cs-", "rayl/re-ff-", new G4LogLogInterpolation, 250*eV, 100*GeV), |
---|
57 | intrinsicLowEnergyLimit(10*eV), |
---|
58 | intrinsicHighEnergyLimit(100*GeV) |
---|
59 | { |
---|
60 | if (GetLowEnergyLimit() < intrinsicLowEnergyLimit || |
---|
61 | GetHighEnergyLimit() > intrinsicHighEnergyLimit) |
---|
62 | G4Exception("G4LowEnergyPolarizedRayleigh::G4LowEnergyPolarizedRayleigh - Energy limit outside intrinsic process validity range"); |
---|
63 | |
---|
64 | |
---|
65 | G4cout << G4endl; |
---|
66 | G4cout << "*******************************************************************************" << G4endl; |
---|
67 | G4cout << "*******************************************************************************" << G4endl; |
---|
68 | G4cout << " The class G4LowEnergyPolarizedRayleigh is NOT SUPPORTED ANYMORE. " << G4endl; |
---|
69 | G4cout << " It will be REMOVED with the next major release of Geant4. " << G4endl; |
---|
70 | G4cout << " Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl; |
---|
71 | G4cout << "*******************************************************************************" << G4endl; |
---|
72 | G4cout << "*******************************************************************************" << G4endl; |
---|
73 | G4cout << G4endl; |
---|
74 | } |
---|
75 | |
---|
76 | |
---|
77 | G4VParticleChange* G4LowEnergyPolarizedRayleigh::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep) |
---|
78 | { |
---|
79 | // aParticleChange comes from G4VProcess |
---|
80 | aParticleChange.Initialize(aTrack); |
---|
81 | |
---|
82 | const G4DynamicParticle* incomingPhoton = aTrack.GetDynamicParticle(); |
---|
83 | G4double incomingPhotonEnergy = incomingPhoton->GetKineticEnergy(); |
---|
84 | |
---|
85 | if (incomingPhotonEnergy <= GetLowEnergyLimit()) |
---|
86 | { |
---|
87 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
88 | aParticleChange.ProposeEnergy(0.); |
---|
89 | aParticleChange.ProposeLocalEnergyDeposit(incomingPhotonEnergy); |
---|
90 | |
---|
91 | return G4VLowEnergyDiscretePhotonProcess::PostStepDoIt(aTrack, aStep); |
---|
92 | } |
---|
93 | |
---|
94 | const G4VCrossSectionHandler* crossSectionHandle = GetCrossSectionHandler(); |
---|
95 | const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple(); |
---|
96 | G4int zAtom = crossSectionHandle->SelectRandomAtom(couple, incomingPhotonEnergy); |
---|
97 | |
---|
98 | G4double outcomingPhotonCosTheta = GenerateCosTheta(incomingPhotonEnergy, zAtom); |
---|
99 | G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta); |
---|
100 | G4double beta=GeneratePolarizationAngle(); |
---|
101 | |
---|
102 | // incomingPhoton reference frame: |
---|
103 | // z = versor parallel to the incomingPhotonDirection |
---|
104 | // x = versor parallel to the incomingPhotonPolarization |
---|
105 | // y = defined as z^x |
---|
106 | |
---|
107 | // outgoingPhoton reference frame: |
---|
108 | // z' = versor parallel to the outgoingPhotonDirection |
---|
109 | // x' = defined as x-x*z'z' normalized |
---|
110 | // y' = defined as z'^x' |
---|
111 | |
---|
112 | G4ThreeVector z(incomingPhoton->GetMomentumDirection().unit()); |
---|
113 | G4ThreeVector x(GetPhotonPolarization(*incomingPhoton)); |
---|
114 | G4ThreeVector y(z.cross(x)); |
---|
115 | |
---|
116 | // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z |
---|
117 | G4double xDir; |
---|
118 | G4double yDir; |
---|
119 | G4double zDir; |
---|
120 | zDir=outcomingPhotonCosTheta; |
---|
121 | xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta); |
---|
122 | yDir=xDir; |
---|
123 | xDir*=std::cos(outcomingPhotonPhi); |
---|
124 | yDir*=std::sin(outcomingPhotonPhi); |
---|
125 | |
---|
126 | G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit()); |
---|
127 | G4ThreeVector xPrime(x.perpPart(zPrime).unit()); |
---|
128 | G4ThreeVector yPrime(zPrime.cross(xPrime)); |
---|
129 | |
---|
130 | // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta) |
---|
131 | G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta)); |
---|
132 | |
---|
133 | aParticleChange.ProposeEnergy(incomingPhotonEnergy); |
---|
134 | aParticleChange.ProposeMomentumDirection(zPrime); |
---|
135 | aParticleChange.ProposePolarization(outcomingPhotonPolarization); |
---|
136 | aParticleChange.SetNumberOfSecondaries(0); |
---|
137 | |
---|
138 | // returns aParticleChange though pParticleChange and G4VProcess::PostStepDoIt |
---|
139 | return G4VLowEnergyDiscretePhotonProcess::PostStepDoIt(aTrack, aStep); |
---|
140 | } |
---|
141 | |
---|
142 | |
---|
143 | |
---|
144 | |
---|
145 | G4double G4LowEnergyPolarizedRayleigh::GenerateCosTheta(G4double incomingPhotonEnergy, G4int zAtom) const |
---|
146 | { |
---|
147 | // d sigma k0 |
---|
148 | // --------- = r0^2 * pi * F^2(x, Z) * ( 2 - sin^2 theta) * std::sin (theta), x = ---- std::sin(theta/2) |
---|
149 | // d theta hc |
---|
150 | |
---|
151 | // d sigma k0 1 - y |
---|
152 | // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + y^2), x = ---- std::sqrt ( ------- ), y = std::cos(theta) |
---|
153 | // d y hc 2 |
---|
154 | |
---|
155 | // Z |
---|
156 | // F(x, Z) ~ -------- |
---|
157 | // a + bx |
---|
158 | // |
---|
159 | // The time to exit from the outer loop grows as ~ k0 |
---|
160 | // On pcgeant2 the time is ~ 1 s for k0 ~ 1 MeV on the oxygen element. A 100 GeV |
---|
161 | // event will take ~ 10 hours. |
---|
162 | // |
---|
163 | // On the avarage the inner loop does 1.5 iterations before exiting |
---|
164 | |
---|
165 | const G4double xFactor = (incomingPhotonEnergy*cm)/(h_Planck*c_light); |
---|
166 | const G4VEMDataSet * formFactorData = GetScatterFunctionData(); |
---|
167 | |
---|
168 | G4double cosTheta; |
---|
169 | G4double fCosTheta; |
---|
170 | G4double x; |
---|
171 | G4double fValue; |
---|
172 | |
---|
173 | do |
---|
174 | { |
---|
175 | do |
---|
176 | { |
---|
177 | cosTheta = 2.*G4UniformRand()-1.; |
---|
178 | fCosTheta = (1.+cosTheta*cosTheta)/2.; |
---|
179 | } |
---|
180 | while (fCosTheta < G4UniformRand()); |
---|
181 | |
---|
182 | x = xFactor*std::sqrt((1.-cosTheta)/2.); |
---|
183 | |
---|
184 | if (x > 1.e+005) |
---|
185 | fValue = formFactorData->FindValue(x, zAtom-1); |
---|
186 | else |
---|
187 | fValue = formFactorData->FindValue(0., zAtom-1); |
---|
188 | |
---|
189 | fValue/=zAtom; |
---|
190 | fValue*=fValue; |
---|
191 | } |
---|
192 | while(fValue < G4UniformRand()); |
---|
193 | |
---|
194 | return cosTheta; |
---|
195 | } |
---|
196 | |
---|
197 | |
---|
198 | |
---|
199 | G4double G4LowEnergyPolarizedRayleigh::GeneratePhi(G4double cosTheta) const |
---|
200 | { |
---|
201 | // d sigma |
---|
202 | // --------- = alpha * ( 1 - sin^2 (theta) * cos^2 (phi) ) |
---|
203 | // d phi |
---|
204 | |
---|
205 | // On the average the loop takes no more than 2 iterations before exiting |
---|
206 | |
---|
207 | G4double phi; |
---|
208 | G4double cosPhi; |
---|
209 | G4double phiProbability; |
---|
210 | G4double sin2Theta; |
---|
211 | |
---|
212 | sin2Theta=1.-cosTheta*cosTheta; |
---|
213 | |
---|
214 | do |
---|
215 | { |
---|
216 | phi = twopi * G4UniformRand(); |
---|
217 | cosPhi = std::cos(phi); |
---|
218 | phiProbability= 1. - sin2Theta*cosPhi*cosPhi; |
---|
219 | } |
---|
220 | while (phiProbability < G4UniformRand()); |
---|
221 | |
---|
222 | return phi; |
---|
223 | } |
---|
224 | |
---|
225 | |
---|
226 | |
---|
227 | |
---|
228 | |
---|
229 | G4double G4LowEnergyPolarizedRayleigh::GeneratePolarizationAngle(void) const |
---|
230 | { |
---|
231 | // Rayleigh polarization is always on the x' direction |
---|
232 | |
---|
233 | return 0; |
---|
234 | } |
---|