source: trunk/source/processes/electromagnetic/lowenergy/src/G4LivermorePolarizedRayleighModel.cc@ 980

Last change on this file since 980 was 968, checked in by garnier, 17 years ago

fichier ajoutes

File size: 12.8 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// $Id: G4LivermorePolarizedRayleighModel.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 "G4LivermorePolarizedRayleighModel.hh"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33
34using namespace std;
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(const G4ParticleDefinition*,
39 const G4String& nam)
40:G4VEmModel(nam),isInitialised(false),crossSectionHandler(0),formFactorData(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 Polarized 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
65G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel()
66{
67 if (crossSectionHandler) delete crossSectionHandler;
68 if (formFactorData) delete formFactorData;
69}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73void G4LivermorePolarizedRayleighModel::Initialise(const G4ParticleDefinition* particle,
74 const G4DataVector& cuts)
75{
76// Rayleigh process: The Quantum Theory of Radiation
77// W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
78// Scattering function: A simple model of photon transport
79// D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
80// Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter
81// T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
82
83 if (verboseLevel > 3)
84 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
85
86 if (crossSectionHandler)
87 {
88 crossSectionHandler->Clear();
89 delete crossSectionHandler;
90 }
91
92 // Energy limits
93
94 if (LowEnergyLimit() < lowEnergyLimit)
95 {
96 G4cout << "G4LivermorePolarizedRayleighModel: low energy limit increased from " <<
97 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
98 SetLowEnergyLimit(lowEnergyLimit);
99 }
100 if (HighEnergyLimit() > highEnergyLimit)
101 {
102 G4cout << "G4LivermorePolarizedRayleighModel: high energy limit decreased from " <<
103 HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
104 SetHighEnergyLimit(highEnergyLimit);
105 }
106
107 // Read data files for all materials
108
109 crossSectionHandler = new G4CrossSectionHandler;
110 crossSectionHandler->Clear();
111 G4String crossSectionFile = "rayl/re-cs-";
112 crossSectionHandler->LoadData(crossSectionFile);
113
114 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
115 G4String formFactorFile = "rayl/re-ff-";
116 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
117 formFactorData->LoadData(formFactorFile);
118
119 //
120 if (verboseLevel > 2)
121 G4cout << "Loaded cross section files for Livermore Polarized Rayleigh model" << G4endl;
122
123 InitialiseElementSelectors(particle,cuts);
124
125 G4cout << "Livermore Polarized Rayleigh model is initialized " << G4endl
126 << "Energy range: "
127 << LowEnergyLimit() / keV << " keV - "
128 << HighEnergyLimit() / GeV << " GeV"
129 << G4endl;
130
131 if(isInitialised) return;
132
133 if(pParticleChange)
134 fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
135 else
136 fParticleChange = new G4ParticleChangeForGamma();
137
138 isInitialised = true;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
143G4double G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(
144 const G4ParticleDefinition*,
145 G4double GammaEnergy,
146 G4double Z, G4double,
147 G4double, G4double)
148{
149 if (verboseLevel > 3)
150 G4cout << "Calling CrossSectionPerAtom() of G4LivermorePolarizedRayleighModel" << G4endl;
151
152 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
153 return cs;
154}
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157
158void G4LivermorePolarizedRayleighModel::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 G4LivermorePolarizedRayleighModel" << 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 // SI - IS THE FOLLOWING RETURN NECESSARY ?
175 return ;
176 }
177
178 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
179
180 // Select randomly one element in the current material
181 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
182
183 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
184 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
185 G4double beta=GeneratePolarizationAngle();
186
187 // incomingPhoton reference frame:
188 // z = versor parallel to the incomingPhotonDirection
189 // x = versor parallel to the incomingPhotonPolarization
190 // y = defined as z^x
191
192 // outgoingPhoton reference frame:
193 // z' = versor parallel to the outgoingPhotonDirection
194 // x' = defined as x-x*z'z' normalized
195 // y' = defined as z'^x'
196
197 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit());
198 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
199 G4ThreeVector y(z.cross(x));
200
201 // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
202 G4double xDir;
203 G4double yDir;
204 G4double zDir;
205 zDir=outcomingPhotonCosTheta;
206 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
207 yDir=xDir;
208 xDir*=std::cos(outcomingPhotonPhi);
209 yDir*=std::sin(outcomingPhotonPhi);
210
211 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
212 G4ThreeVector xPrime(x.perpPart(zPrime).unit());
213 G4ThreeVector yPrime(zPrime.cross(xPrime));
214
215 // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta)
216 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
217
218 fParticleChange->ProposeMomentumDirection(zPrime);
219 fParticleChange->ProposePolarization(outcomingPhotonPolarization);
220 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
221
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
225
226G4double G4LivermorePolarizedRayleighModel::GenerateCosTheta(G4double incomingPhotonEnergy, G4int zAtom) const
227{
228 // d sigma k0
229 // --------- = r0^2 * pi * F^2(x, Z) * ( 2 - sin^2 theta) * std::sin (theta), x = ---- std::sin(theta/2)
230 // d theta hc
231
232 // d sigma k0 1 - y
233 // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + y^2), x = ---- std::sqrt ( ------- ), y = std::cos(theta)
234 // d y hc 2
235
236 // Z
237 // F(x, Z) ~ --------
238 // a + bx
239 //
240 // The time to exit from the outer loop grows as ~ k0
241 // On pcgeant2 the time is ~ 1 s for k0 ~ 1 MeV on the oxygen element. A 100 GeV
242 // event will take ~ 10 hours.
243 //
244 // On the avarage the inner loop does 1.5 iterations before exiting
245
246 const G4double xFactor = (incomingPhotonEnergy*cm)/(h_Planck*c_light);
247 //const G4VEMDataSet * formFactorData = GetScatterFunctionData();
248
249 G4double cosTheta;
250 G4double fCosTheta;
251 G4double x;
252 G4double fValue;
253
254 do
255 {
256 do
257 {
258 cosTheta = 2.*G4UniformRand()-1.;
259 fCosTheta = (1.+cosTheta*cosTheta)/2.;
260 }
261 while (fCosTheta < G4UniformRand());
262
263 x = xFactor*std::sqrt((1.-cosTheta)/2.);
264
265 if (x > 1.e+005)
266 fValue = formFactorData->FindValue(x, zAtom-1);
267 else
268 fValue = formFactorData->FindValue(0., zAtom-1);
269
270 fValue/=zAtom;
271 fValue*=fValue;
272 }
273 while(fValue < G4UniformRand());
274
275 return cosTheta;
276}
277
278//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279
280G4double G4LivermorePolarizedRayleighModel::GeneratePhi(G4double cosTheta) const
281{
282 // d sigma
283 // --------- = alpha * ( 1 - sin^2 (theta) * cos^2 (phi) )
284 // d phi
285
286 // On the average the loop takes no more than 2 iterations before exiting
287
288 G4double phi;
289 G4double cosPhi;
290 G4double phiProbability;
291 G4double sin2Theta;
292
293 sin2Theta=1.-cosTheta*cosTheta;
294
295 do
296 {
297 phi = twopi * G4UniformRand();
298 cosPhi = std::cos(phi);
299 phiProbability= 1. - sin2Theta*cosPhi*cosPhi;
300 }
301 while (phiProbability < G4UniformRand());
302
303 return phi;
304}
305
306//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
307
308G4double G4LivermorePolarizedRayleighModel::GeneratePolarizationAngle(void) const
309{
310 // Rayleigh polarization is always on the x' direction
311
312 return 0;
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316
317G4ThreeVector G4LivermorePolarizedRayleighModel::GetPhotonPolarization(const G4DynamicParticle& photon)
318{
319
320// SI - From G4VLowEnergyDiscretePhotonProcess.cc
321
322 G4ThreeVector photonMomentumDirection;
323 G4ThreeVector photonPolarization;
324
325 photonPolarization = photon.GetPolarization();
326 photonMomentumDirection = photon.GetMomentumDirection();
327
328 if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
329 {
330 // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
331 // then polarization is choosen randomly.
332
333 G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
334 G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
335
336 G4double angle(G4UniformRand() * twopi);
337
338 e1*=std::cos(angle);
339 e2*=std::sin(angle);
340
341 photonPolarization=e1+e2;
342 }
343 else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
344 {
345 // if |photonPolarization * photonDirection0| != 0.
346 // then polarization is made orthonormal;
347
348 photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
349 }
350
351 return photonPolarization.unit();
352}
353
Note: See TracBrowser for help on using the repository browser.