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