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 // ------------------------------------------- 26 // -------------------------------------------------------------------- 27 // 27 // 28 // 28 // 29 // Author: A. Forti 29 // Author: A. Forti 30 // Maria Grazia Pia (Maria.Grazia.Pia@ 30 // Maria Grazia Pia (Maria.Grazia.Pia@cern.ch) 31 // 31 // 32 // History: 32 // History: 33 // -------- 33 // -------- 34 // Added Livermore data table construction met 34 // Added Livermore data table construction methods A. Forti 35 // Added BuildMeanFreePath A. Forti 35 // Added BuildMeanFreePath A. Forti 36 // Added PostStepDoIt A. Forti 36 // Added PostStepDoIt A. Forti 37 // Added SelectRandomAtom A. Forti 37 // Added SelectRandomAtom A. Forti 38 // Added map of the elements A.Forti 38 // Added map of the elements A.Forti 39 // 24.04.01 V.Ivanchenko remove RogueWave 39 // 24.04.01 V.Ivanchenko remove RogueWave 40 // 11.08.2001 MGP - Major revision according t 40 // 11.08.2001 MGP - Major revision according to a design iteration 41 // 06.10.2001 MGP - Added strategy to test ran 41 // 06.10.2001 MGP - Added strategy to test range for secondary generation 42 // 05.06.2002 F.Longo and G.Depaola - bug fix 42 // 05.06.2002 F.Longo and G.Depaola - bug fixed in angular distribution 43 // 20.10.2002 G. Depaola - Change sampling m 43 // 20.10.2002 G. Depaola - Change sampling method of theta 44 // 22.01.2003 V.Ivanchenko - Cut per region 44 // 22.01.2003 V.Ivanchenko - Cut per region 45 // 24.04.2003 V.Ivanchenko - Cut per region mf 45 // 24.04.2003 V.Ivanchenko - Cut per region mfpt 46 // 46 // 47 // ------------------------------------------- 47 // -------------------------------------------------------------------- 48 48 49 #include "G4LowEnergyRayleigh.hh" 49 #include "G4LowEnergyRayleigh.hh" 50 #include "Randomize.hh" 50 #include "Randomize.hh" 51 #include "G4PhysicalConstants.hh" 51 #include "G4PhysicalConstants.hh" 52 #include "G4SystemOfUnits.hh" 52 #include "G4SystemOfUnits.hh" 53 #include "G4ParticleDefinition.hh" 53 #include "G4ParticleDefinition.hh" 54 #include "G4Track.hh" 54 #include "G4Track.hh" 55 #include "G4Step.hh" 55 #include "G4Step.hh" 56 #include "G4ForceCondition.hh" 56 #include "G4ForceCondition.hh" 57 #include "G4Gamma.hh" 57 #include "G4Gamma.hh" 58 #include "G4Electron.hh" 58 #include "G4Electron.hh" 59 #include "G4DynamicParticle.hh" 59 #include "G4DynamicParticle.hh" 60 #include "G4VParticleChange.hh" 60 #include "G4VParticleChange.hh" 61 #include "G4ThreeVector.hh" 61 #include "G4ThreeVector.hh" 62 #include "G4RDVCrossSectionHandler.hh" 62 #include "G4RDVCrossSectionHandler.hh" 63 #include "G4RDCrossSectionHandler.hh" 63 #include "G4RDCrossSectionHandler.hh" 64 #include "G4RDVEMDataSet.hh" 64 #include "G4RDVEMDataSet.hh" 65 #include "G4RDCompositeEMDataSet.hh" 65 #include "G4RDCompositeEMDataSet.hh" 66 #include "G4RDVDataSetAlgorithm.hh" 66 #include "G4RDVDataSetAlgorithm.hh" 67 #include "G4RDLogLogInterpolation.hh" 67 #include "G4RDLogLogInterpolation.hh" 68 68 69 #include "G4MaterialCutsCouple.hh" 69 #include "G4MaterialCutsCouple.hh" 70 70 71 G4LowEnergyRayleigh::G4LowEnergyRayleigh(const 71 G4LowEnergyRayleigh::G4LowEnergyRayleigh(const G4String& processName) 72 : G4VDiscreteProcess(processName), 72 : G4VDiscreteProcess(processName), 73 lowEnergyLimit(250*eV), 73 lowEnergyLimit(250*eV), 74 highEnergyLimit(100*GeV), 74 highEnergyLimit(100*GeV), 75 intrinsicLowEnergyLimit(10*eV), 75 intrinsicLowEnergyLimit(10*eV), 76 intrinsicHighEnergyLimit(100*GeV) 76 intrinsicHighEnergyLimit(100*GeV) 77 { 77 { 78 if (lowEnergyLimit < intrinsicLowEnergyLimit 78 if (lowEnergyLimit < intrinsicLowEnergyLimit || 79 highEnergyLimit > intrinsicHighEnergyLim 79 highEnergyLimit > intrinsicHighEnergyLimit) 80 { 80 { 81 G4Exception("G4LowEnergyRayleigh::G4LowE 81 G4Exception("G4LowEnergyRayleigh::G4LowEnergyRayleigh()", 82 "OutOfRange", FatalException 82 "OutOfRange", FatalException, 83 "Energy limit outside intrin 83 "Energy limit outside intrinsic process validity range!"); 84 } 84 } 85 85 86 crossSectionHandler = new G4RDCrossSectionHa 86 crossSectionHandler = new G4RDCrossSectionHandler(); 87 87 88 G4RDVDataSetAlgorithm* ffInterpolation = new 88 G4RDVDataSetAlgorithm* ffInterpolation = new G4RDLogLogInterpolation; 89 G4String formFactorFile = "rayl/re-ff-"; 89 G4String formFactorFile = "rayl/re-ff-"; 90 formFactorData = new G4RDCompositeEMDataSet( 90 formFactorData = new G4RDCompositeEMDataSet(ffInterpolation,1.,1.); 91 formFactorData->LoadData(formFactorFile); 91 formFactorData->LoadData(formFactorFile); 92 92 93 meanFreePathTable = 0; 93 meanFreePathTable = 0; 94 94 95 if (verboseLevel > 0) 95 if (verboseLevel > 0) 96 { 96 { 97 G4cout << GetProcessName() << " is crea 97 G4cout << GetProcessName() << " is created " << G4endl 98 << "Energy range: " 98 << "Energy range: " 99 << lowEnergyLimit / keV << " keV - " 99 << lowEnergyLimit / keV << " keV - " 100 << highEnergyLimit / GeV << " GeV" 100 << highEnergyLimit / GeV << " GeV" 101 << G4endl; 101 << G4endl; 102 } 102 } 103 } 103 } 104 104 105 G4LowEnergyRayleigh::~G4LowEnergyRayleigh() 105 G4LowEnergyRayleigh::~G4LowEnergyRayleigh() 106 { 106 { 107 delete meanFreePathTable; 107 delete meanFreePathTable; 108 delete crossSectionHandler; 108 delete crossSectionHandler; 109 delete formFactorData; 109 delete formFactorData; 110 } 110 } 111 111 112 void G4LowEnergyRayleigh::BuildPhysicsTable(co 112 void G4LowEnergyRayleigh::BuildPhysicsTable(const G4ParticleDefinition& ) 113 { 113 { 114 114 115 crossSectionHandler->Clear(); 115 crossSectionHandler->Clear(); 116 G4String crossSectionFile = "rayl/re-cs-"; 116 G4String crossSectionFile = "rayl/re-cs-"; 117 crossSectionHandler->LoadData(crossSectionFi 117 crossSectionHandler->LoadData(crossSectionFile); 118 118 119 delete meanFreePathTable; 119 delete meanFreePathTable; 120 meanFreePathTable = crossSectionHandler->Bui 120 meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials(); 121 } 121 } 122 122 123 G4VParticleChange* G4LowEnergyRayleigh::PostSt 123 G4VParticleChange* G4LowEnergyRayleigh::PostStepDoIt(const G4Track& aTrack, 124 const G4Step& aStep) 124 const G4Step& aStep) 125 { 125 { 126 126 127 aParticleChange.Initialize(aTrack); 127 aParticleChange.Initialize(aTrack); 128 128 129 const G4DynamicParticle* incidentPhoton = aT 129 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle(); 130 G4double photonEnergy0 = incidentPhoton->Get 130 G4double photonEnergy0 = incidentPhoton->GetKineticEnergy(); 131 131 132 if (photonEnergy0 <= lowEnergyLimit) 132 if (photonEnergy0 <= lowEnergyLimit) 133 { 133 { 134 aParticleChange.ProposeTrackStatus(fStop 134 aParticleChange.ProposeTrackStatus(fStopAndKill); 135 aParticleChange.ProposeEnergy(0.); 135 aParticleChange.ProposeEnergy(0.); 136 aParticleChange.ProposeLocalEnergyDeposi 136 aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0); 137 return G4VDiscreteProcess::PostStepDoIt( 137 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep); 138 } 138 } 139 139 140 // G4double e0m = photonEnergy0 / electron_ 140 // G4double e0m = photonEnergy0 / electron_mass_c2 ; 141 G4ParticleMomentum photonDirection0 = incide 141 G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection(); 142 142 143 // Select randomly one element in the curren 143 // Select randomly one element in the current material 144 const G4MaterialCutsCouple* couple = aTrack. 144 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple(); 145 G4int Z = crossSectionHandler->SelectRandomA 145 G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0); 146 146 147 // Sample the angle of the scattered photon 147 // Sample the angle of the scattered photon 148 148 149 G4double wlPhoton = h_Planck*c_light/photonE 149 G4double wlPhoton = h_Planck*c_light/photonEnergy0; 150 150 151 G4double gReject,x,dataFormFactor; 151 G4double gReject,x,dataFormFactor; 152 G4double randomFormFactor; 152 G4double randomFormFactor; 153 G4double cosTheta; 153 G4double cosTheta; 154 G4double sinTheta; 154 G4double sinTheta; 155 G4double fcostheta; 155 G4double fcostheta; 156 156 157 do 157 do 158 { 158 { 159 do 159 do 160 { 160 { 161 cosTheta = 2. * G4UniformRand() - 1.; 161 cosTheta = 2. * G4UniformRand() - 1.; 162 fcostheta = ( 1. + cosTheta*cosTheta)/2. 162 fcostheta = ( 1. + cosTheta*cosTheta)/2.; 163 } while (fcostheta < G4UniformRand()); 163 } while (fcostheta < G4UniformRand()); 164 164 165 G4double sinThetaHalf = std::sqrt((1. - 165 G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.); 166 x = sinThetaHalf / (wlPhoton/cm); 166 x = sinThetaHalf / (wlPhoton/cm); 167 if (x > 1.e+005) 167 if (x > 1.e+005) 168 dataFormFactor = formFactorData->Find 168 dataFormFactor = formFactorData->FindValue(x,Z-1); 169 else 169 else 170 dataFormFactor = formFactorData->Find 170 dataFormFactor = formFactorData->FindValue(0.,Z-1); 171 randomFormFactor = G4UniformRand() * Z * 171 randomFormFactor = G4UniformRand() * Z * Z; 172 sinTheta = std::sqrt(1. - cosTheta*cosTh 172 sinTheta = std::sqrt(1. - cosTheta*cosTheta); 173 gReject = dataFormFactor * dataFormFacto 173 gReject = dataFormFactor * dataFormFactor; 174 174 175 } while( gReject < randomFormFactor); 175 } while( gReject < randomFormFactor); 176 176 177 // Scattered photon angles. ( Z - axis along 177 // Scattered photon angles. ( Z - axis along the parent photon) 178 G4double phi = twopi * G4UniformRand() ; 178 G4double phi = twopi * G4UniformRand() ; 179 G4double dirX = sinTheta*std::cos(phi); 179 G4double dirX = sinTheta*std::cos(phi); 180 G4double dirY = sinTheta*std::sin(phi); 180 G4double dirY = sinTheta*std::sin(phi); 181 G4double dirZ = cosTheta; 181 G4double dirZ = cosTheta; 182 182 183 // Update G4VParticleChange for the scattere 183 // Update G4VParticleChange for the scattered photon 184 G4ThreeVector photonDirection1(dirX, dirY, d 184 G4ThreeVector photonDirection1(dirX, dirY, dirZ); 185 185 186 photonDirection1.rotateUz(photonDirection0); 186 photonDirection1.rotateUz(photonDirection0); 187 aParticleChange.ProposeEnergy(photonEnergy0) 187 aParticleChange.ProposeEnergy(photonEnergy0); 188 aParticleChange.ProposeMomentumDirection(pho 188 aParticleChange.ProposeMomentumDirection(photonDirection1); 189 189 190 aParticleChange.SetNumberOfSecondaries(0); 190 aParticleChange.SetNumberOfSecondaries(0); 191 191 192 return G4VDiscreteProcess::PostStepDoIt(aTra 192 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 193 } 193 } 194 194 195 G4bool G4LowEnergyRayleigh::IsApplicable(const 195 G4bool G4LowEnergyRayleigh::IsApplicable(const G4ParticleDefinition& particle) 196 { 196 { 197 return ( &particle == G4Gamma::Gamma() ); 197 return ( &particle == G4Gamma::Gamma() ); 198 } 198 } 199 199 200 G4double G4LowEnergyRayleigh::GetMeanFreePath( 200 G4double G4LowEnergyRayleigh::GetMeanFreePath(const G4Track& track, 201 G4double, // previousStepSize 201 G4double, // previousStepSize 202 G4ForceCondition*) 202 G4ForceCondition*) 203 { 203 { 204 const G4DynamicParticle* photon = track.GetD 204 const G4DynamicParticle* photon = track.GetDynamicParticle(); 205 G4double energy = photon->GetKineticEnergy() 205 G4double energy = photon->GetKineticEnergy(); 206 const G4MaterialCutsCouple* couple = track.G 206 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple(); 207 size_t materialIndex = couple->GetIndex(); 207 size_t materialIndex = couple->GetIndex(); 208 208 209 G4double meanFreePath; 209 G4double meanFreePath; 210 if (energy > highEnergyLimit) meanFreePath = 210 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex); 211 else if (energy < lowEnergyLimit) meanFreePa 211 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX; 212 else meanFreePath = meanFreePathTable->FindV 212 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex); 213 return meanFreePath; 213 return meanFreePath; 214 } 214 } 215 215