| 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.5 2009/05/02 15:20:53 sincerti Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-03-cand-01 $
|
|---|
| 28 | //
|
|---|
| 29 | // History:
|
|---|
| 30 | // --------
|
|---|
| 31 | // 02 May 2009 S Incerti as V. Ivanchenko proposed in G4LivermoreRayleighModel.cc
|
|---|
| 32 | //
|
|---|
| 33 | // Cleanup initialisation and generation of secondaries:
|
|---|
| 34 | // - apply internal high-energy limit only in constructor
|
|---|
| 35 | // - do not apply low-energy limit (default is 0)
|
|---|
| 36 | // - remove GetMeanFreePath method and table
|
|---|
| 37 | // - remove initialisation of element selector
|
|---|
| 38 | // - use G4ElementSelector
|
|---|
| 39 |
|
|---|
| 40 | #include "G4LivermorePolarizedRayleighModel.hh"
|
|---|
| 41 |
|
|---|
| 42 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 43 |
|
|---|
| 44 | using namespace std;
|
|---|
| 45 |
|
|---|
| 46 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 47 |
|
|---|
| 48 | G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(const G4ParticleDefinition*,
|
|---|
| 49 | const G4String& nam)
|
|---|
| 50 | :G4VEmModel(nam),isInitialised(false),crossSectionHandler(0),formFactorData(0)
|
|---|
| 51 | {
|
|---|
| 52 | lowEnergyLimit = 250 * eV;
|
|---|
| 53 | highEnergyLimit = 100 * GeV;
|
|---|
| 54 |
|
|---|
| 55 | //SetLowEnergyLimit(lowEnergyLimit);
|
|---|
| 56 | SetHighEnergyLimit(highEnergyLimit);
|
|---|
| 57 | //
|
|---|
| 58 | verboseLevel= 0;
|
|---|
| 59 | // Verbosity scale:
|
|---|
| 60 | // 0 = nothing
|
|---|
| 61 | // 1 = warning for energy non-conservation
|
|---|
| 62 | // 2 = details of energy budget
|
|---|
| 63 | // 3 = calculation of cross sections, file openings, sampling of atoms
|
|---|
| 64 | // 4 = entering in methods
|
|---|
| 65 |
|
|---|
| 66 | if(verboseLevel > 0) {
|
|---|
| 67 | G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl
|
|---|
| 68 | << "Energy range: "
|
|---|
| 69 | << lowEnergyLimit / eV << " eV - "
|
|---|
| 70 | << highEnergyLimit / GeV << " GeV"
|
|---|
| 71 | << G4endl;
|
|---|
| 72 | }
|
|---|
| 73 | }
|
|---|
| 74 |
|
|---|
| 75 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 76 |
|
|---|
| 77 | G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel()
|
|---|
| 78 | {
|
|---|
| 79 | if (crossSectionHandler) delete crossSectionHandler;
|
|---|
| 80 | if (formFactorData) delete formFactorData;
|
|---|
| 81 | }
|
|---|
| 82 |
|
|---|
| 83 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 84 |
|
|---|
| 85 | void G4LivermorePolarizedRayleighModel::Initialise(const G4ParticleDefinition* particle,
|
|---|
| 86 | const G4DataVector& cuts)
|
|---|
| 87 | {
|
|---|
| 88 | // Rayleigh process: The Quantum Theory of Radiation
|
|---|
| 89 | // W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
|
|---|
| 90 | // Scattering function: A simple model of photon transport
|
|---|
| 91 | // D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
|
|---|
| 92 | // Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter
|
|---|
| 93 | // T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
|
|---|
| 94 |
|
|---|
| 95 | if (verboseLevel > 3)
|
|---|
| 96 | G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
|
|---|
| 97 |
|
|---|
| 98 | if (crossSectionHandler)
|
|---|
| 99 | {
|
|---|
| 100 | crossSectionHandler->Clear();
|
|---|
| 101 | delete crossSectionHandler;
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 | // Read data files for all materials
|
|---|
| 105 |
|
|---|
| 106 | crossSectionHandler = new G4CrossSectionHandler;
|
|---|
| 107 | crossSectionHandler->Clear();
|
|---|
| 108 | G4String crossSectionFile = "rayl/re-cs-";
|
|---|
| 109 | crossSectionHandler->LoadData(crossSectionFile);
|
|---|
| 110 |
|
|---|
| 111 | G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
|
|---|
| 112 | G4String formFactorFile = "rayl/re-ff-";
|
|---|
| 113 | formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
|
|---|
| 114 | formFactorData->LoadData(formFactorFile);
|
|---|
| 115 |
|
|---|
| 116 | InitialiseElementSelectors(particle,cuts);
|
|---|
| 117 |
|
|---|
| 118 | //
|
|---|
| 119 | if (verboseLevel > 2)
|
|---|
| 120 | G4cout << "Loaded cross section files for Livermore Polarized Rayleigh model" << G4endl;
|
|---|
| 121 |
|
|---|
| 122 | InitialiseElementSelectors(particle,cuts);
|
|---|
| 123 |
|
|---|
| 124 | if (verboseLevel > 0) {
|
|---|
| 125 | G4cout << "Livermore Polarized Rayleigh model is initialized " << G4endl
|
|---|
| 126 | << "Energy range: "
|
|---|
| 127 | << LowEnergyLimit() / eV << " eV - "
|
|---|
| 128 | << HighEnergyLimit() / GeV << " GeV"
|
|---|
| 129 | << G4endl;
|
|---|
| 130 | }
|
|---|
| 131 |
|
|---|
| 132 | if(isInitialised) return;
|
|---|
| 133 | fParticleChange = GetParticleChangeForGamma();
|
|---|
| 134 | isInitialised = true;
|
|---|
| 135 | }
|
|---|
| 136 |
|
|---|
| 137 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 138 |
|
|---|
| 139 | G4double G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(
|
|---|
| 140 | const G4ParticleDefinition*,
|
|---|
| 141 | G4double GammaEnergy,
|
|---|
| 142 | G4double Z, G4double,
|
|---|
| 143 | G4double, G4double)
|
|---|
| 144 | {
|
|---|
| 145 | if (verboseLevel > 3)
|
|---|
| 146 | G4cout << "Calling CrossSectionPerAtom() of G4LivermorePolarizedRayleighModel" << G4endl;
|
|---|
| 147 |
|
|---|
| 148 | if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0.0;
|
|---|
| 149 |
|
|---|
| 150 | G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
|
|---|
| 151 | return cs;
|
|---|
| 152 | }
|
|---|
| 153 |
|
|---|
| 154 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 155 |
|
|---|
| 156 | void G4LivermorePolarizedRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
|
|---|
| 157 | const G4MaterialCutsCouple* couple,
|
|---|
| 158 | const G4DynamicParticle* aDynamicGamma,
|
|---|
| 159 | G4double,
|
|---|
| 160 | G4double)
|
|---|
| 161 | {
|
|---|
| 162 | if (verboseLevel > 3)
|
|---|
| 163 | G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
|
|---|
| 164 |
|
|---|
| 165 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
|
|---|
| 166 |
|
|---|
| 167 | if (photonEnergy0 <= lowEnergyLimit)
|
|---|
| 168 | {
|
|---|
| 169 | fParticleChange->ProposeTrackStatus(fStopAndKill);
|
|---|
| 170 | fParticleChange->SetProposedKineticEnergy(0.);
|
|---|
| 171 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
|
|---|
| 172 | return ;
|
|---|
| 173 | }
|
|---|
| 174 |
|
|---|
| 175 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
|
|---|
| 176 |
|
|---|
| 177 | // Select randomly one element in the current material
|
|---|
| 178 | // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
|
|---|
| 179 | const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
|
|---|
| 180 | const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
|
|---|
| 181 | G4int Z = (G4int)elm->GetZ();
|
|---|
| 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 |
|
|---|
| 226 | G4double 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 |
|
|---|
| 280 | G4double 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 |
|
|---|
| 308 | G4double 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 |
|
|---|
| 317 | G4ThreeVector 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 |
|
|---|