Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // >> 26 // $Id: G4LivermorePolarizedRayleighModel.cc,v 1.1 2008/10/30 14:16:35 sincerti Exp $ >> 27 // GEANT4 tag $Name: geant4-09-02-patch-03 $ 26 // 28 // 27 // Author: Sebastien Incerti << 28 // 30 October 2008 << 29 // on base of G4LowEnergyPolarizedRayl << 30 // << 31 // History: << 32 // -------- << 33 // 02 May 2009 S Incerti as V. Ivanchenko pr << 34 // << 35 // Cleanup initialisation and generation of se << 36 // - apply internal high-ener << 37 // - do not apply low-energy << 38 // - remove GetMeanFreePath m << 39 // - remove initialisation of << 40 // - use G4ElementSelector << 41 29 42 #include "G4LivermorePolarizedRayleighModel.hh 30 #include "G4LivermorePolarizedRayleighModel.hh" 43 #include "G4PhysicalConstants.hh" << 44 #include "G4SystemOfUnits.hh" << 45 #include "G4LogLogInterpolation.hh" << 46 #include "G4CompositeEMDataSet.hh" << 47 #include "G4AutoLock.hh" << 48 31 49 //....oooOO0OOooo........oooOO0OOooo........oo 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 50 33 51 using namespace std; 34 using namespace std; 52 namespace { G4Mutex LivermorePolarizedRayleigh << 53 35 54 //....oooOO0OOooo........oooOO0OOooo........oo 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 55 37 56 G4PhysicsFreeVector* G4LivermorePolarizedRayle << 57 G4PhysicsFreeVector* G4LivermorePolarizedRayle << 58 << 59 G4LivermorePolarizedRayleighModel::G4Livermore 38 G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(const G4ParticleDefinition*, 60 const G4String& nam) << 39 const G4String& nam) 61 :G4VEmModel(nam),fParticleChange(nullptr),is << 40 :G4VEmModel(nam),isInitialised(false) 62 { 41 { 63 lowEnergyLimit = 250 * CLHEP::eV; << 42 lowEnergyLimit = 250 * eV; // SI - Could be 10 eV ? >> 43 highEnergyLimit = 100 * GeV; >> 44 >> 45 SetLowEnergyLimit(lowEnergyLimit); >> 46 SetHighEnergyLimit(highEnergyLimit); 64 // 47 // 65 verboseLevel= 0; 48 verboseLevel= 0; 66 // Verbosity scale: 49 // Verbosity scale: 67 // 0 = nothing 50 // 0 = nothing 68 // 1 = warning for energy non-conservation 51 // 1 = warning for energy non-conservation 69 // 2 = details of energy budget 52 // 2 = details of energy budget 70 // 3 = calculation of cross sections, file o 53 // 3 = calculation of cross sections, file openings, sampling of atoms 71 // 4 = entering in methods 54 // 4 = entering in methods 72 55 73 if(verboseLevel > 0) { << 56 G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl 74 G4cout << "Livermore Polarized Rayleigh is << 57 << "Energy range: " 75 << "Energy range: " << LowEnergyLim << 58 << lowEnergyLimit / keV << " keV - " 76 << HighEnergyLimit() / CLHEP::GeV << " Ge << 59 << highEnergyLimit / GeV << " GeV" 77 << G4endl; << 60 << G4endl; 78 } << 79 } 61 } 80 62 81 //....oooOO0OOooo........oooOO0OOooo........oo 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 82 64 83 G4LivermorePolarizedRayleighModel::~G4Livermor 65 G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel() 84 { 66 { 85 if(IsMaster()) { << 67 delete crossSectionHandler; 86 for(G4int i=0; i<maxZ; ++i) { << 68 delete formFactorData; 87 if(dataCS[i]) { << 88 delete dataCS[i]; << 89 dataCS[i] = nullptr; << 90 delete formFactorData[i]; << 91 formFactorData[i] = nullptr; << 92 } << 93 } << 94 } << 95 } 69 } 96 70 97 //....oooOO0OOooo........oooOO0OOooo........oo 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 98 72 99 void G4LivermorePolarizedRayleighModel::Initia 73 void G4LivermorePolarizedRayleighModel::Initialise(const G4ParticleDefinition* particle, 100 const G 74 const G4DataVector& cuts) 101 { 75 { 102 // Rayleigh process: Th << 76 // Rayleigh process: The Quantum Theory of Radiation 103 // W. << 77 // W. Heitler, Oxford at the Clarendon Press, Oxford (1954) 104 // Scattering function: A << 78 // Scattering function: A simple model of photon transport 105 // D. << 79 // D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510 106 // Polarization of the outcoming photon: Be << 80 // Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter 107 // X- << 81 // T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168 108 // T. << 82 109 if (verboseLevel > 3) 83 if (verboseLevel > 3) 110 G4cout << "Calling G4LivermorePolarizedRay 84 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl; 111 85 112 if(IsMaster()) { << 86 InitialiseElementSelectors(particle,cuts); 113 << 114 // Initialise element selector << 115 InitialiseElementSelectors(particle, cuts) << 116 << 117 // Access to elements << 118 const char* path = G4FindDataDir("G4LEDATA << 119 auto elmTable = G4Element::GetElementTable << 120 for (auto const & elm : *elmTable) { << 121 G4int Z = std::min(elm->GetZasInt(), max << 122 if( nullptr == dataCS[Z] ) { ReadData(Z, << 123 } << 124 } << 125 << 126 if(isInitialised) { return; } << 127 fParticleChange = GetParticleChangeForGamma( << 128 isInitialised = true; << 129 } << 130 << 131 << 132 //....oooOO0OOooo........oooOO0OOooo........oo << 133 << 134 void G4LivermorePolarizedRayleighModel::Initia << 135 const G4ParticleDefinition*, G << 136 { << 137 SetElementSelectors(masterModel->GetElementS << 138 } << 139 87 140 //....oooOO0OOooo........oooOO0OOooo........oo << 88 // Energy limits 141 << 142 void G4LivermorePolarizedRayleighModel::ReadDa << 143 { << 144 if (verboseLevel > 1) { << 145 G4cout << "Calling ReadData() of G4Livermo << 146 << G4endl; << 147 } << 148 << 149 if(nullptr != dataCS[Z]) { return; } << 150 << 151 const char* datadir = path; << 152 89 153 if(nullptr == datadir) << 90 if (LowEnergyLimit() < lowEnergyLimit) 154 { 91 { 155 datadir = G4FindDataDir("G4LEDATA"); << 92 G4cout << "G4LivermorePolarizedRayleighModel: low energy limit increased from " << 156 if(nullptr == datadir) << 93 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl; 157 { << 94 SetLowEnergyLimit(lowEnergyLimit); 158 G4Exception("G4LivermoreRayleighModelModel << 159 FatalException, << 160 "Environment variable G4LEDATA not d << 161 return; << 162 } << 163 } 95 } 164 dataCS[Z] = new G4PhysicsFreeVector(); << 96 if (HighEnergyLimit() > highEnergyLimit) 165 formFactorData[Z] = new G4PhysicsFreeVector( << 97 { 166 << 98 G4cout << "G4LivermorePolarizedRayleighModel: high energy limit decreased from " << 167 std::ostringstream ostCS; << 99 HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl; 168 ostCS << datadir << "/livermore/rayl/re-cs-" << 100 SetHighEnergyLimit(highEnergyLimit); 169 std::ifstream finCS(ostCS.str().c_str()); << 170 << 171 if( !finCS .is_open() ) { << 172 G4ExceptionDescription ed; << 173 ed << "G4LivermorePolarizedRayleighModel d << 174 << "> is not opened!" << G4endl; << 175 G4Exception("G4LivermorePolarizedRayleighM << 176 FatalException, << 177 ed,"G4LEDATA version should be G4EMLOW8.0 << 178 return; << 179 } else { << 180 if(verboseLevel > 3) { << 181 G4cout << "File " << ostCS.str() << 182 << " is opened by G4LivermoreRayleighMo << 183 } 101 } 184 dataCS[Z]->Retrieve(finCS, true); << 102 185 } << 103 // Read data files for all materials 186 104 187 std::ostringstream ostFF; << 105 crossSectionHandler = new G4CrossSectionHandler; 188 ostFF << datadir << "/livermore/rayl/re-ff-" << 106 crossSectionHandler->Clear(); 189 std::ifstream finFF(ostFF.str().c_str()); << 107 G4String crossSectionFile = "rayl/re-cs-"; 190 << 108 crossSectionHandler->LoadData(crossSectionFile); 191 if( !finFF.is_open() ) { << 109 192 G4ExceptionDescription ed; << 110 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation; 193 ed << "G4LivermorePolarizedRayleighModel d << 111 G4String formFactorFile = "rayl/re-ff-"; 194 << "> is not opened!" << G4endl; << 112 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.); 195 G4Exception("G4LivermorePolarizedRayleighM << 113 formFactorData->LoadData(formFactorFile); 196 FatalException, << 114 197 ed,"G4LEDATA version should be G4EMLOW8.0 << 115 // 198 return; << 116 if (verboseLevel > 2) 199 } else { << 117 G4cout << "Loaded cross section files for Livermore Polarized Rayleigh model" << G4endl; 200 if(verboseLevel > 3) { << 118 201 G4cout << "File " << ostFF.str() << 119 G4cout << "Livermore Polarized Rayleigh model is initialized " << G4endl 202 << " is opened by G4LivermoreRa << 120 << "Energy range: " 203 } << 121 << LowEnergyLimit() / keV << " keV - " 204 formFactorData[Z]->Retrieve(finFF, true); << 122 << HighEnergyLimit() / GeV << " GeV" 205 } << 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; 206 } 133 } 207 << 134 208 //....oooOO0OOooo........oooOO0OOooo........oo 135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 209 << 136 210 G4double G4LivermorePolarizedRayleighModel::Co 137 G4double G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom( 211 const << 138 const G4ParticleDefinition*, 212 G4double GammaEnergy, << 139 G4double GammaEnergy, 213 G4double Z, G4double, << 140 G4double Z, G4double, 214 G4double, G4double) << 141 G4double, G4double) 215 { 142 { 216 if (verboseLevel > 1) { << 143 if (verboseLevel > 3) 217 G4cout << "G4LivermoreRayleighModel::Compu << 144 G4cout << "Calling CrossSectionPerAtom() of G4LivermorePolarizedRayleighModel" << G4endl; 218 << G4endl; << 145 219 } << 146 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy); 220 << 147 return cs; 221 if(GammaEnergy < lowEnergyLimit) { return 0. << 222 << 223 G4double xs = 0.0; << 224 << 225 G4int intZ = G4lrint(Z); << 226 if(intZ < 1 || intZ > maxZ) { return xs; } << 227 << 228 G4PhysicsFreeVector* pv = dataCS[intZ]; << 229 << 230 // if element was not initialised << 231 // do initialisation safely for MT mode << 232 if(nullptr == pv) { << 233 InitialiseForElement(0, intZ); << 234 pv = dataCS[intZ]; << 235 if(nullptr == pv) { return xs; } << 236 } << 237 << 238 G4int n = G4int(pv->GetVectorLength() - 1); << 239 G4double e = GammaEnergy/MeV; << 240 if(e >= pv->Energy(n)) { << 241 xs = (*pv)[n]/(e*e); << 242 } else if(e >= pv->Energy(0)) { << 243 xs = pv->Value(e)/(e*e); << 244 } << 245 return xs; << 246 } 148 } 247 149 248 //....oooOO0OOooo........oooOO0OOooo........oo 150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 249 151 250 void G4LivermorePolarizedRayleighModel::Sample << 152 void G4LivermorePolarizedRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/, 251 std::vector<G4DynamicParticle* << 153 const G4MaterialCutsCouple* couple, 252 const G4MaterialCutsCouple* couple, << 154 const G4DynamicParticle* aDynamicGamma, 253 const G4DynamicParticle* aDynamicGamma, << 155 G4double, 254 G4double, G4double) << 156 G4double) 255 { 157 { 256 if (verboseLevel > 3) 158 if (verboseLevel > 3) 257 G4cout << "Calling SampleSecondaries() of 159 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl; 258 160 259 G4double photonEnergy0 = aDynamicGamma->GetK 161 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); 260 162 261 if (photonEnergy0 <= lowEnergyLimit) 163 if (photonEnergy0 <= lowEnergyLimit) 262 { 164 { 263 fParticleChange->ProposeTrackStatus(fStopA << 165 fParticleChange->ProposeTrackStatus(fStopAndKill); 264 fParticleChange->SetProposedKineticEnergy( << 166 fParticleChange->SetProposedKineticEnergy(0.); 265 fParticleChange->ProposeLocalEnergyDeposit << 167 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); 266 return; << 168 // SI - IS THE FOLLOWING RETURN NECESSARY ? >> 169 return ; 267 } 170 } 268 171 269 G4ParticleMomentum photonDirection0 = aDynam 172 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); 270 173 271 // Select randomly one element in the curren 174 // Select randomly one element in the current material 272 const G4ParticleDefinition* particle = aDyn << 175 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0); 273 const G4Element* elm = SelectRandomAtom(coup << 274 G4int Z = elm->GetZasInt(); << 275 176 276 G4double outcomingPhotonCosTheta = GenerateC 177 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z); 277 G4double outcomingPhotonPhi = GeneratePhi(ou 178 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta); 278 G4double beta = GeneratePolarizationAngle(); << 179 G4double beta=GeneratePolarizationAngle(); 279 180 280 // incomingPhoton reference frame: 181 // incomingPhoton reference frame: 281 // z = versor parallel to the incomingPhoton 182 // z = versor parallel to the incomingPhotonDirection 282 // x = versor parallel to the incomingPhoton 183 // x = versor parallel to the incomingPhotonPolarization 283 // y = defined as z^x 184 // y = defined as z^x 284 185 285 // outgoingPhoton reference frame: 186 // outgoingPhoton reference frame: 286 // z' = versor parallel to the outgoingPhoto 187 // z' = versor parallel to the outgoingPhotonDirection 287 // x' = defined as x-x*z'z' normalized 188 // x' = defined as x-x*z'z' normalized 288 // y' = defined as z'^x' << 189 // y' = defined as z'^x' >> 190 289 G4ThreeVector z(aDynamicGamma->GetMomentumDi 191 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit()); 290 G4ThreeVector x(GetPhotonPolarization(*aDyna 192 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma)); 291 G4ThreeVector y(z.cross(x)); 193 G4ThreeVector y(z.cross(x)); 292 194 293 // z' = std::cos(phi)*std::sin(theta) << 195 // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z 294 // x + std::sin(phi)*std::sin(theta) y + std << 295 G4double xDir; 196 G4double xDir; 296 G4double yDir; 197 G4double yDir; 297 G4double zDir; 198 G4double zDir; 298 zDir=outcomingPhotonCosTheta; 199 zDir=outcomingPhotonCosTheta; 299 xDir=std::sqrt(1-outcomingPhotonCosTheta*out 200 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta); 300 yDir=xDir; 201 yDir=xDir; 301 xDir*=std::cos(outcomingPhotonPhi); 202 xDir*=std::cos(outcomingPhotonPhi); 302 yDir*=std::sin(outcomingPhotonPhi); 203 yDir*=std::sin(outcomingPhotonPhi); 303 204 304 G4ThreeVector zPrime((xDir*x + yDir*y + zDir 205 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit()); 305 G4ThreeVector xPrime(x.perpPart(zPrime).unit 206 G4ThreeVector xPrime(x.perpPart(zPrime).unit()); 306 G4ThreeVector yPrime(zPrime.cross(xPrime)); 207 G4ThreeVector yPrime(zPrime.cross(xPrime)); 307 208 308 // outgoingPhotonPolarization is directed as << 209 // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta) 309 // x' std::cos(beta) + y' std::sin(beta) << 310 G4ThreeVector outcomingPhotonPolarization(xP 210 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta)); 311 211 312 fParticleChange->ProposeMomentumDirection(zP 212 fParticleChange->ProposeMomentumDirection(zPrime); 313 fParticleChange->ProposePolarization(outcomi 213 fParticleChange->ProposePolarization(outcomingPhotonPolarization); 314 fParticleChange->SetProposedKineticEnergy(ph 214 fParticleChange->SetProposedKineticEnergy(photonEnergy0); >> 215 315 } 216 } 316 217 317 //....oooOO0OOooo........oooOO0OOooo........oo 218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 318 219 319 G4double G4LivermorePolarizedRayleighModel::Ge << 220 G4double G4LivermorePolarizedRayleighModel::GenerateCosTheta(G4double incomingPhotonEnergy, G4int zAtom) const 320 { 221 { 321 // d sigma 222 // d sigma k0 322 // --------- = r0^2 * pi * F^2(x, Z) * ( 2 223 // --------- = r0^2 * pi * F^2(x, Z) * ( 2 - sin^2 theta) * std::sin (theta), x = ---- std::sin(theta/2) 323 // d theta 224 // d theta hc 324 225 325 // d sigma 226 // d sigma k0 1 - y 326 // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + 227 // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + y^2), x = ---- std::sqrt ( ------- ), y = std::cos(theta) 327 // d y 228 // d y hc 2 328 229 329 // Z 230 // Z 330 // F(x, Z) ~ -------- 231 // F(x, Z) ~ -------- 331 // a + bx 232 // a + bx 332 // 233 // 333 // The time to exit from the outer loop grow 234 // The time to exit from the outer loop grows as ~ k0 334 // On pcgeant2 the time is ~ 1 s for k0 ~ 1 235 // On pcgeant2 the time is ~ 1 s for k0 ~ 1 MeV on the oxygen element. A 100 GeV 335 // event will take ~ 10 hours. 236 // event will take ~ 10 hours. 336 // 237 // 337 // On the avarage the inner loop does 1.5 it 238 // On the avarage the inner loop does 1.5 iterations before exiting 338 const G4double xxfact = CLHEP::cm/(CLHEP::h_ << 239 339 const G4double xFactor = incomingPhotonEnerg << 240 const G4double xFactor = (incomingPhotonEnergy*cm)/(h_Planck*c_light); >> 241 //const G4VEMDataSet * formFactorData = GetScatterFunctionData(); 340 242 341 G4double cosTheta; 243 G4double cosTheta; 342 G4double fCosTheta; 244 G4double fCosTheta; 343 G4double x; 245 G4double x; 344 G4double fValue; 246 G4double fValue; 345 247 346 if (incomingPhotonEnergy > 5.*CLHEP::MeV) << 248 do 347 { << 348 cosTheta = 1.; << 349 } << 350 else << 351 { << 352 do << 353 { 249 { 354 do 250 do 355 { 251 { 356 cosTheta = 2.*G4UniformRand()-1.; 252 cosTheta = 2.*G4UniformRand()-1.; 357 fCosTheta = (1.+cosTheta*cosTheta)/2.; 253 fCosTheta = (1.+cosTheta*cosTheta)/2.; 358 } 254 } 359 while (fCosTheta < G4UniformRand()); 255 while (fCosTheta < G4UniformRand()); 360 256 361 x = xFactor*std::sqrt((1.-cosTheta)/2.); 257 x = xFactor*std::sqrt((1.-cosTheta)/2.); 362 258 363 if (x > 1.e+005) 259 if (x > 1.e+005) 364 fValue = formFactorData[Z]->Value(x); << 260 fValue = formFactorData->FindValue(x, zAtom-1); 365 else 261 else 366 fValue = formFactorData[Z]->Value(0.); << 262 fValue = formFactorData->FindValue(0., zAtom-1); 367 263 368 fValue /= Z; << 264 fValue/=zAtom; 369 fValue *= fValue; << 265 fValue*=fValue; 370 } 266 } 371 while(fValue < G4UniformRand()); << 267 while(fValue < G4UniformRand()); 372 } << 373 268 374 return cosTheta; 269 return cosTheta; 375 } 270 } 376 271 377 //....oooOO0OOooo........oooOO0OOooo........oo 272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 378 273 379 G4double G4LivermorePolarizedRayleighModel::Ge 274 G4double G4LivermorePolarizedRayleighModel::GeneratePhi(G4double cosTheta) const 380 { 275 { 381 // d sigma 276 // d sigma 382 // --------- = alpha * ( 1 - sin^2 (theta) * 277 // --------- = alpha * ( 1 - sin^2 (theta) * cos^2 (phi) ) 383 // d phi 278 // d phi 384 279 385 // On the average the loop takes no more tha 280 // On the average the loop takes no more than 2 iterations before exiting 386 281 387 G4double phi; 282 G4double phi; 388 G4double cosPhi; 283 G4double cosPhi; 389 G4double phiProbability; 284 G4double phiProbability; 390 G4double sin2Theta; 285 G4double sin2Theta; 391 286 392 sin2Theta=1.-cosTheta*cosTheta; 287 sin2Theta=1.-cosTheta*cosTheta; 393 288 394 do 289 do 395 { 290 { 396 phi = CLHEP::twopi * G4UniformRand(); << 291 phi = twopi * G4UniformRand(); 397 cosPhi = std::cos(phi); 292 cosPhi = std::cos(phi); 398 phiProbability= 1. - sin2Theta*cosPhi*co 293 phiProbability= 1. - sin2Theta*cosPhi*cosPhi; 399 } 294 } 400 while (phiProbability < G4UniformRand()); 295 while (phiProbability < G4UniformRand()); 401 296 402 return phi; 297 return phi; 403 } 298 } 404 299 405 //....oooOO0OOooo........oooOO0OOooo........oo 300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 406 301 407 G4double G4LivermorePolarizedRayleighModel::Ge 302 G4double G4LivermorePolarizedRayleighModel::GeneratePolarizationAngle(void) const 408 { 303 { 409 // Rayleigh polarization is always on the x' 304 // Rayleigh polarization is always on the x' direction >> 305 410 return 0; 306 return 0; 411 } 307 } 412 308 413 //....oooOO0OOooo........oooOO0OOooo........oo 309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 414 310 415 G4ThreeVector G4LivermorePolarizedRayleighMode 311 G4ThreeVector G4LivermorePolarizedRayleighModel::GetPhotonPolarization(const G4DynamicParticle& photon) 416 { 312 { 417 // From G4VLowEnergyDiscretePhotonProcess.cc << 313 >> 314 // SI - From G4VLowEnergyDiscretePhotonProcess.cc >> 315 418 G4ThreeVector photonMomentumDirection; 316 G4ThreeVector photonMomentumDirection; 419 G4ThreeVector photonPolarization; 317 G4ThreeVector photonPolarization; 420 318 421 photonPolarization = photon.GetPolarization( 319 photonPolarization = photon.GetPolarization(); 422 photonMomentumDirection = photon.GetMomentum 320 photonMomentumDirection = photon.GetMomentumDirection(); 423 321 424 if ((!photonPolarization.isOrthogonal(photon 322 if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.) 425 { 323 { 426 // if |photonPolarization|==0. or |photo 324 // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0| 427 // then polarization is choosen randomly << 325 // then polarization is choosen randomly. >> 326 428 G4ThreeVector e1(photonMomentumDirection 327 G4ThreeVector e1(photonMomentumDirection.orthogonal().unit()); 429 G4ThreeVector e2(photonMomentumDirection 328 G4ThreeVector e2(photonMomentumDirection.cross(e1).unit()); 430 329 431 G4double angle(G4UniformRand() * CLHEP:: << 330 G4double angle(G4UniformRand() * twopi); 432 331 433 e1*=std::cos(angle); 332 e1*=std::cos(angle); 434 e2*=std::sin(angle); 333 e2*=std::sin(angle); 435 334 436 photonPolarization=e1+e2; 335 photonPolarization=e1+e2; 437 } 336 } 438 else if (photonPolarization.howOrthogonal(ph 337 else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.) 439 { 338 { 440 // if |photonPolarization * photonDirect 339 // if |photonPolarization * photonDirection0| != 0. 441 // then polarization is made orthonormal << 340 // then polarization is made orthonormal; >> 341 442 photonPolarization=photonPolarization.pe 342 photonPolarization=photonPolarization.perpPart(photonMomentumDirection); 443 } 343 } 444 344 445 return photonPolarization.unit(); 345 return photonPolarization.unit(); 446 } 346 } 447 347 448 //....oooOO0OOooo........oooOO0OOooo........oo << 449 << 450 void G4LivermorePolarizedRayleighModel::Initia << 451 const G4ParticleDefinition*, G << 452 { << 453 G4AutoLock l(&LivermorePolarizedRayleighMode << 454 if(nullptr == dataCS[Z]) { ReadData(Z); } << 455 l.unlock(); << 456 } << 457 348