[819] | 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 | // -------------------------------------------------------------------- |
---|
| 27 | // |
---|
[1192] | 28 | // $Id: G4PenelopeRayleigh.cc,v 1.19 2009/06/11 15:47:08 mantero Exp $ |
---|
[1228] | 29 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[819] | 30 | // |
---|
| 31 | // Author: L. Pandola (luciano.pandola@cern.ch) |
---|
| 32 | // |
---|
| 33 | // History: |
---|
| 34 | // -------- |
---|
| 35 | // 14 Feb 2003 MG Pia Corrected compilation errors and warnings |
---|
| 36 | // from SUN |
---|
| 37 | // 10 Mar 2003 V.Ivanchenko Remove CutPerMaterial warning |
---|
| 38 | // 12 Mar 2003 L.Pandola Code "cleaned" - Cuts per region |
---|
| 39 | // 17 Mar 2004 L.Pandola Removed unnecessary calls to std::pow(a,b) |
---|
| 40 | // 18 Mar 2004 M.Mendenhall Introduced SamplingTable (performance improvement) |
---|
| 41 | // 03 Sep 2007 L.Pandola Bug fix for the filling of physics table for |
---|
| 42 | // compounds defined by the mass fraction (bug #965) |
---|
[1055] | 43 | // 02 Apr 2009 L.Pandola Bux fixed in the calculation of mfp for compound |
---|
| 44 | // materials defined as fraction of mass |
---|
| 45 | // (reported by Zhang Qiwei) |
---|
[819] | 46 | // -------------------------------------------------------------------- |
---|
| 47 | |
---|
| 48 | #include "G4PenelopeRayleigh.hh" |
---|
| 49 | #include "Randomize.hh" |
---|
| 50 | #include "G4ParticleDefinition.hh" |
---|
| 51 | #include "G4Track.hh" |
---|
| 52 | #include "G4Step.hh" |
---|
| 53 | #include "G4ForceCondition.hh" |
---|
| 54 | #include "G4Gamma.hh" |
---|
| 55 | #include "G4Electron.hh" |
---|
| 56 | #include "G4DynamicParticle.hh" |
---|
| 57 | #include "G4VParticleChange.hh" |
---|
| 58 | #include "G4ThreeVector.hh" |
---|
| 59 | #include "G4VCrossSectionHandler.hh" |
---|
| 60 | #include "G4CrossSectionHandler.hh" |
---|
| 61 | #include "G4VEMDataSet.hh" |
---|
| 62 | #include "G4EMDataSet.hh" |
---|
| 63 | #include "G4CompositeEMDataSet.hh" |
---|
| 64 | #include "G4VDataSetAlgorithm.hh" |
---|
| 65 | #include "G4LogLogInterpolation.hh" |
---|
| 66 | #include "G4PenelopeIntegrator.hh" |
---|
| 67 | #include "G4MaterialCutsCouple.hh" |
---|
| 68 | |
---|
| 69 | G4PenelopeRayleigh::G4PenelopeRayleigh(const G4String& processName) |
---|
| 70 | : G4VDiscreteProcess(processName), |
---|
| 71 | lowEnergyLimit(250*eV), |
---|
| 72 | highEnergyLimit(100*GeV), |
---|
| 73 | samplingConstant(0.0), |
---|
| 74 | nBins(200), |
---|
| 75 | intrinsicLowEnergyLimit(10*eV), |
---|
| 76 | intrinsicHighEnergyLimit(100*GeV) |
---|
| 77 | |
---|
| 78 | { |
---|
| 79 | if (lowEnergyLimit < intrinsicLowEnergyLimit || |
---|
| 80 | highEnergyLimit > intrinsicHighEnergyLimit) |
---|
| 81 | { |
---|
| 82 | G4Exception("G4PenelopeRayleigh::G4PenelopeRayleigh - energy limit outside intrinsic process validity range"); |
---|
| 83 | } |
---|
| 84 | |
---|
| 85 | samplingFunction_x = 0; |
---|
| 86 | samplingFunction_y = 0; |
---|
| 87 | meanFreePathTable = 0; |
---|
| 88 | material = 0; |
---|
| 89 | |
---|
| 90 | if (verboseLevel > 0) |
---|
| 91 | { |
---|
| 92 | G4cout << GetProcessName() << " is created " << G4endl |
---|
| 93 | << "Energy range: " |
---|
| 94 | << lowEnergyLimit / keV << " keV - " |
---|
| 95 | << highEnergyLimit / GeV << " GeV" |
---|
| 96 | << G4endl; |
---|
| 97 | } |
---|
[1055] | 98 | |
---|
| 99 | G4cout << G4endl; |
---|
| 100 | G4cout << "*******************************************************************************" << G4endl; |
---|
| 101 | G4cout << "*******************************************************************************" << G4endl; |
---|
| 102 | G4cout << " The class G4PenelopeRayleigh is NOT SUPPORTED ANYMORE. " << G4endl; |
---|
| 103 | G4cout << " It will be REMOVED with the next major release of Geant4. " << G4endl; |
---|
| 104 | G4cout << " Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl; |
---|
| 105 | G4cout << "*******************************************************************************" << G4endl; |
---|
| 106 | G4cout << "*******************************************************************************" << G4endl; |
---|
| 107 | G4cout << G4endl; |
---|
[819] | 108 | } |
---|
| 109 | |
---|
| 110 | G4PenelopeRayleigh::~G4PenelopeRayleigh() |
---|
| 111 | { |
---|
| 112 | delete meanFreePathTable; |
---|
| 113 | |
---|
| 114 | SamplingTablePair::iterator i; |
---|
| 115 | for(i=SamplingTables.begin(); i != SamplingTables.end(); i++) { |
---|
| 116 | delete (*i).second.first; |
---|
| 117 | delete (*i).second.second; |
---|
| 118 | } |
---|
| 119 | } |
---|
| 120 | |
---|
| 121 | |
---|
| 122 | void G4PenelopeRayleigh::BuildPhysicsTable(const G4ParticleDefinition& ) |
---|
| 123 | { |
---|
| 124 | |
---|
| 125 | G4DataVector energyVector; |
---|
| 126 | G4double dBin = std::log10(highEnergyLimit/lowEnergyLimit)/nBins; |
---|
| 127 | for (G4int i=0;i<nBins;i++) |
---|
| 128 | { |
---|
| 129 | energyVector.push_back(std::pow(10.,std::log10(lowEnergyLimit)+i*dBin)); |
---|
| 130 | } |
---|
| 131 | |
---|
| 132 | const G4MaterialTable* materialTable = G4Material::GetMaterialTable(); |
---|
| 133 | G4int nMaterials = G4Material::GetNumberOfMaterials(); |
---|
| 134 | |
---|
| 135 | size_t nOfBins = energyVector.size(); |
---|
| 136 | size_t bin=0; |
---|
| 137 | |
---|
| 138 | G4VDataSetAlgorithm* algo = new G4LogLogInterpolation(); |
---|
| 139 | G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo,1.,1.); |
---|
| 140 | std::vector<G4VEMDataSet*> matCrossSections; |
---|
| 141 | |
---|
| 142 | G4int m; |
---|
| 143 | for (m=0; m<nMaterials; m++) |
---|
| 144 | { |
---|
| 145 | G4DataVector* energies = new G4DataVector; |
---|
| 146 | G4DataVector* data = new G4DataVector; |
---|
| 147 | material = (*materialTable)[m]; |
---|
| 148 | |
---|
| 149 | G4int nElements = material->GetNumberOfElements(); |
---|
| 150 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
| 151 | const G4double k1=849.3315; |
---|
| 152 | |
---|
| 153 | |
---|
| 154 | G4int IZZ=0; |
---|
| 155 | G4int iright=0; |
---|
| 156 | for (G4int i=0; i<nElements; i++) { |
---|
| 157 | G4int Z = (G4int) (*elementVector)[i]->GetZ(); |
---|
| 158 | if (Z>IZZ){ |
---|
| 159 | IZZ = Z; |
---|
| 160 | iright=i; |
---|
| 161 | } |
---|
| 162 | } |
---|
| 163 | for (bin=0; bin<nOfBins; bin++) |
---|
| 164 | { |
---|
| 165 | energies->push_back(energyVector[bin]); |
---|
| 166 | G4double ec=std::min(energyVector[bin],0.5*IZZ); |
---|
| 167 | G4double energyRatio = ec/electron_mass_c2; |
---|
| 168 | facte = k1*energyRatio*energyRatio; |
---|
| 169 | G4double cs=0; |
---|
| 170 | G4PenelopeIntegrator<G4PenelopeRayleigh,G4double(G4PenelopeRayleigh::*)(G4double)> theIntegrator; |
---|
| 171 | cs = |
---|
| 172 | theIntegrator.Calculate(this,&G4PenelopeRayleigh::DifferentialCrossSection,-1.0,0.90,1e-06); |
---|
| 173 | cs += theIntegrator.Calculate(this,&G4PenelopeRayleigh::DifferentialCrossSection,0.90,0.9999999,1e-06); |
---|
| 174 | cs = cs*(ec/energyVector[bin])*(ec/energyVector[bin])*pi*classic_electr_radius*classic_electr_radius; |
---|
| 175 | const G4double* vector_of_atoms = material->GetVecNbOfAtomsPerVolume(); |
---|
| 176 | const G4int* stechiometric = material->GetAtomsVector(); |
---|
[1055] | 177 | //cs is the cross section _per atom_ in the case of compounds, while it is |
---|
| 178 | //_per molecule_ in the case of molecules |
---|
| 179 | G4double density = material->GetTotNbOfAtomsPerVolume(); //non-bound molecules (default) |
---|
[819] | 180 | if (stechiometric) |
---|
| 181 | { |
---|
| 182 | if (stechiometric[iright]) |
---|
| 183 | density = vector_of_atoms[iright]/stechiometric[iright]; //number of molecules per volume |
---|
| 184 | } |
---|
| 185 | G4double cross = density*cs; |
---|
| 186 | data->push_back(cross); |
---|
| 187 | } |
---|
| 188 | G4VEMDataSet* elSet = new G4EMDataSet(0,energies,data,algo); |
---|
| 189 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo); |
---|
| 190 | setForMat->AddComponent(elSet); |
---|
| 191 | matCrossSections.push_back(setForMat); |
---|
| 192 | } |
---|
| 193 | G4double matCS = 0.0; |
---|
| 194 | |
---|
| 195 | for (m=0; m<nMaterials; m++) |
---|
| 196 | { |
---|
| 197 | G4DataVector* energies = new G4DataVector; |
---|
| 198 | G4DataVector* data = new G4DataVector; |
---|
| 199 | for (bin=0;bin<nOfBins;bin++){ |
---|
| 200 | energies->push_back(energyVector[bin]); |
---|
| 201 | matCS = (matCrossSections[m]->GetComponent(0))->FindValue(energyVector[bin]); |
---|
| 202 | if (matCS > 0.){ |
---|
| 203 | data->push_back(1./matCS); |
---|
| 204 | } |
---|
| 205 | else |
---|
| 206 | { |
---|
| 207 | data->push_back(DBL_MAX); |
---|
| 208 | } |
---|
| 209 | } |
---|
| 210 | G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.); |
---|
| 211 | materialSet->AddComponent(dataSet); |
---|
| 212 | } |
---|
| 213 | meanFreePathTable = materialSet; |
---|
| 214 | } |
---|
| 215 | |
---|
| 216 | G4VParticleChange* G4PenelopeRayleigh::PostStepDoIt(const G4Track& aTrack, |
---|
| 217 | const G4Step& aStep) |
---|
| 218 | { |
---|
| 219 | aParticleChange.Initialize(aTrack); |
---|
| 220 | |
---|
| 221 | const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle(); |
---|
| 222 | G4double photonEnergy0 = incidentPhoton->GetKineticEnergy(); |
---|
| 223 | |
---|
| 224 | if (photonEnergy0 <= lowEnergyLimit) |
---|
| 225 | { |
---|
| 226 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
| 227 | aParticleChange.ProposeEnergy(0.); |
---|
| 228 | aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0); |
---|
| 229 | return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep); |
---|
| 230 | } |
---|
| 231 | |
---|
| 232 | |
---|
| 233 | G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection(); |
---|
| 234 | const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple(); |
---|
| 235 | material = couple->GetMaterial(); |
---|
| 236 | |
---|
| 237 | // Sampling inizialitation (build internal table) |
---|
| 238 | |
---|
| 239 | InizialiseSampling(); |
---|
| 240 | |
---|
| 241 | // Sample the angle of the scattered photon |
---|
| 242 | const G4double xpar=41.2148; |
---|
| 243 | G4double x2max = 2.0*std::log(xpar*photonEnergy0/electron_mass_c2); |
---|
| 244 | G4int jm; |
---|
| 245 | G4int asize = samplingFunction_x->size(); |
---|
| 246 | if (x2max<(*samplingFunction_x)[1]) |
---|
| 247 | { |
---|
| 248 | jm=0; |
---|
| 249 | } |
---|
| 250 | else if(x2max>(*samplingFunction_x)[asize-2]) |
---|
| 251 | { |
---|
| 252 | jm=asize-2; |
---|
| 253 | } |
---|
| 254 | else |
---|
| 255 | { |
---|
| 256 | jm=(G4int) ((x2max-(*samplingFunction_x)[0])/samplingConstant); |
---|
| 257 | } |
---|
| 258 | G4double rumax = (*samplingFunction_y)[jm]+((*samplingFunction_y)[jm+1]-(*samplingFunction_y)[jm])* |
---|
| 259 | (x2max-(*samplingFunction_x)[jm])/((*samplingFunction_x)[jm+1]-(*samplingFunction_x)[jm]); |
---|
| 260 | G4int j,ju,jt; |
---|
| 261 | G4double ru,denomin,x2rat; |
---|
| 262 | G4double CDT,G,rand; |
---|
| 263 | do{ |
---|
| 264 | ru = rumax + std::log(G4UniformRand()); |
---|
| 265 | j=0; |
---|
| 266 | ju=jm+1; |
---|
| 267 | do{ |
---|
| 268 | jt=(j+ju)/2; //bipartition |
---|
| 269 | if (ru > (*samplingFunction_y)[jt]) |
---|
| 270 | { |
---|
| 271 | j=jt; |
---|
| 272 | } |
---|
| 273 | else |
---|
| 274 | { |
---|
| 275 | ju=jt; |
---|
| 276 | } |
---|
| 277 | }while ((ju-j)>1); |
---|
| 278 | denomin = (*samplingFunction_y)[j+1]-(*samplingFunction_y)[j]; |
---|
| 279 | if (denomin > 1e-12) |
---|
| 280 | { |
---|
| 281 | x2rat = (*samplingFunction_x)[j]+(((*samplingFunction_x)[j+1]-(*samplingFunction_x)[j])* |
---|
| 282 | (ru-(*samplingFunction_y)[j])/denomin)-x2max; |
---|
| 283 | } |
---|
| 284 | else |
---|
| 285 | { |
---|
| 286 | x2rat = (*samplingFunction_x)[j]-x2max; |
---|
| 287 | } |
---|
| 288 | CDT = 1.0-2.0*std::exp(x2rat); |
---|
| 289 | G = 0.5*(1.0+CDT*CDT); |
---|
| 290 | rand = G4UniformRand(); |
---|
| 291 | }while (rand>G); |
---|
| 292 | |
---|
| 293 | G4double cosTheta = CDT; |
---|
| 294 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
| 295 | |
---|
| 296 | |
---|
| 297 | |
---|
| 298 | |
---|
| 299 | // Scattered photon angles. ( Z - axis along the parent photon) |
---|
| 300 | G4double phi = twopi * G4UniformRand() ; |
---|
| 301 | G4double dirX = sinTheta*std::cos(phi); |
---|
| 302 | G4double dirY = sinTheta*std::sin(phi); |
---|
| 303 | G4double dirZ = cosTheta; |
---|
| 304 | |
---|
| 305 | // Update G4VParticleChange for the scattered photon |
---|
| 306 | G4ThreeVector photonDirection1(dirX, dirY, dirZ); |
---|
| 307 | |
---|
| 308 | photonDirection1.rotateUz(photonDirection0); |
---|
| 309 | aParticleChange.ProposeEnergy(photonEnergy0); |
---|
| 310 | aParticleChange.ProposeMomentumDirection(photonDirection1); |
---|
| 311 | |
---|
| 312 | aParticleChange.SetNumberOfSecondaries(0); |
---|
| 313 | |
---|
| 314 | return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); |
---|
| 315 | } |
---|
| 316 | |
---|
| 317 | G4bool G4PenelopeRayleigh::IsApplicable(const G4ParticleDefinition& particle) |
---|
| 318 | { |
---|
| 319 | return ( &particle == G4Gamma::Gamma() ); |
---|
| 320 | } |
---|
| 321 | |
---|
| 322 | G4double G4PenelopeRayleigh::GetMeanFreePath(const G4Track& track, |
---|
| 323 | G4double, // previousStepSize |
---|
| 324 | G4ForceCondition*) |
---|
| 325 | { |
---|
| 326 | const G4DynamicParticle* photon = track.GetDynamicParticle(); |
---|
| 327 | G4double energy = photon->GetKineticEnergy(); |
---|
| 328 | material = track.GetMaterial(); |
---|
| 329 | size_t materialIndex = material->GetIndex(); |
---|
| 330 | |
---|
| 331 | G4double meanFreePath; |
---|
| 332 | if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex); |
---|
| 333 | else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX; |
---|
| 334 | else |
---|
| 335 | { |
---|
| 336 | meanFreePath = meanFreePathTable->FindValue(energy,materialIndex); |
---|
| 337 | } |
---|
| 338 | return meanFreePath; |
---|
| 339 | } |
---|
| 340 | |
---|
| 341 | void G4PenelopeRayleigh::InizialiseSampling() |
---|
| 342 | { |
---|
| 343 | SamplingTablePair::iterator theTable=SamplingTables.find(material); |
---|
| 344 | const G4int points=241; |
---|
| 345 | G4double Xlow = 0.; |
---|
| 346 | G4double Xhigh=1e-04; |
---|
| 347 | G4double fact = std::pow((1e06/Xhigh),(1/240.0)); |
---|
| 348 | samplingConstant=std::log(fact); |
---|
| 349 | |
---|
| 350 | if (theTable==SamplingTables.end()) { //material not inizialized yet |
---|
| 351 | samplingFunction_x = new G4DataVector(); |
---|
| 352 | samplingFunction_y = new G4DataVector(); |
---|
| 353 | G4double sum = 0.0; |
---|
| 354 | G4PenelopeIntegrator<G4PenelopeRayleigh,G4double(G4PenelopeRayleigh::*)(G4double)> theIntegrator; |
---|
| 355 | sum = theIntegrator.Calculate(this,&G4PenelopeRayleigh::MolecularFormFactor, |
---|
| 356 | Xlow,Xhigh,1e-10); |
---|
| 357 | samplingFunction_x->push_back(Xhigh); |
---|
| 358 | samplingFunction_y->push_back(sum); |
---|
| 359 | G4int i; |
---|
| 360 | for (i=1;i<points;i++){ |
---|
| 361 | Xlow=Xhigh; |
---|
| 362 | Xhigh=Xhigh*fact; |
---|
| 363 | sum = theIntegrator.Calculate(this, |
---|
| 364 | &G4PenelopeRayleigh::MolecularFormFactor, |
---|
| 365 | Xlow,Xhigh,1e-10); |
---|
| 366 | samplingFunction_x->push_back(Xhigh); |
---|
| 367 | samplingFunction_y->push_back(sum+(*samplingFunction_y)[i-1]); |
---|
| 368 | } |
---|
| 369 | for (i=0;i<points;i++){ |
---|
| 370 | (*samplingFunction_x)[i]=std::log((*samplingFunction_x)[i]); |
---|
| 371 | (*samplingFunction_y)[i]=std::log((*samplingFunction_y)[i]); |
---|
| 372 | } |
---|
| 373 | SamplingTables[material] = std::pair<G4DataVector*,G4DataVector*> (samplingFunction_x,samplingFunction_y); |
---|
| 374 | } |
---|
| 375 | else { //material already inizialized |
---|
| 376 | samplingFunction_x=(*theTable).second.first; |
---|
| 377 | samplingFunction_y=(*theTable).second.second; |
---|
| 378 | } |
---|
| 379 | } |
---|
| 380 | |
---|
| 381 | |
---|
| 382 | G4double G4PenelopeRayleigh::MolecularFormFactor(G4double y) |
---|
| 383 | { |
---|
| 384 | const G4int ntot=95; |
---|
| 385 | G4double RA1[ntot] = {0.0e0, 3.9265e+0, 4.3100e+1, 5.2757e+1, 2.5021e+1, |
---|
| 386 | 1.2211e+1, 9.3229e+0, 3.2455e+0, 2.4197e+0, 1.5985e+0, |
---|
| 387 | 3.0926e+1, 1.5315e+1, 7.7061e+0, 3.9493e+0, 2.2042e+0, |
---|
| 388 | 1.9453e+1, 1.9354e+1, 8.0374e+0, 8.3779e+1, 5.7370e+1, |
---|
| 389 | 5.2310e+1, 4.7514e+1, 4.3785e+1, 4.2048e+1, 3.6972e+1, |
---|
| 390 | 3.3849e+1, 3.1609e+1, 2.8763e+1, 2.7217e+1, 2.4263e+1, |
---|
| 391 | 2.2403e+1, 1.8606e+1, 1.5143e+1, 1.4226e+1, 1.1792e+1, |
---|
| 392 | 9.7574e+0, 1.2796e+1, 1.2854e+1, 1.2368e+1, 1.0208e+1, |
---|
| 393 | 8.2823e+0, 7.4677e+0, 7.6028e+0, 6.1090e+0, 5.5346e+0, |
---|
| 394 | 4.2340e+0, 4.0444e+0, 4.2905e+0, 4.7950e+0, 5.1112e+0, |
---|
| 395 | 5.2407e+0, 5.2153e+0, 5.1639e+0, 4.8814e+0, 5.8054e+0, |
---|
| 396 | 6.6724e+0, 6.5104e+0, 6.3364e+0, 6.2889e+0, 6.3028e+0, |
---|
| 397 | 6.3853e+0, 6.3475e+0, 6.5779e+0, 6.8486e+0, 7.0993e+0, |
---|
| 398 | 7.6122e+0, 7.9681e+0, 8.3481e+0, 6.3875e+0, 8.0042e+0, |
---|
| 399 | 8.0820e+0, 7.6940e+0, 7.1927e+0, 6.6751e+0, 6.1623e+0, |
---|
| 400 | 5.8335e+0, 5.5599e+0, 4.6551e+0, 4.4327e+0, 4.7601e+0, |
---|
| 401 | 5.2872e+0, 5.6084e+0, 5.7680e+0, 5.8041e+0, 5.7566e+0, |
---|
| 402 | 5.6541e+0, 6.3932e+0, 6.9313e+0, 7.0027e+0, 6.8796e+0, |
---|
| 403 | 6.4739e+0, 6.2405e+0, 6.0081e+0, 5.5708e+0, 5.3680e+0}; |
---|
| 404 | |
---|
| 405 | |
---|
| 406 | G4double RA2[ntot] = {0.0e0, 1.3426e-1, 9.4875e+1,-1.0896e+2,-4.5494e+1, |
---|
| 407 | -1.9572e+1,-1.2382e+1,-3.6827e+0,-2.4542e+0,-1.4453e+0, |
---|
| 408 | 1.3401e+2, 7.9717e+1, 6.2164e+1, 4.0300e+1, 3.1682e+1, |
---|
| 409 | -1.3639e+1,-1.5950e+1,-5.1523e+0, 1.8351e+2, 1.2205e+2, |
---|
| 410 | 1.0007e+2, 8.5632e+1, 7.9145e+1, 6.3675e+1, 6.2954e+1, |
---|
| 411 | 5.6601e+1, 5.4171e+1, 4.8752e+1, 3.8062e+1, 3.9933e+1, |
---|
| 412 | 4.8343e+1, 4.2137e+1, 3.4617e+1, 2.9430e+1, 2.4010e+1, |
---|
| 413 | 1.9744e+1, 4.0009e+1, 5.1614e+1, 5.0456e+1, 3.9088e+1, |
---|
| 414 | 2.6824e+1, 2.2953e+1, 2.4773e+1, 1.6893e+1, 1.4548e+1, |
---|
| 415 | 9.7226e+0, 1.0192e+1, 1.1153e+1, 1.3188e+1, 1.4733e+1, |
---|
| 416 | 1.5644e+1, 1.5939e+1, 1.5923e+1, 1.5254e+1, 2.0748e+1, |
---|
| 417 | 2.6901e+1, 2.7032e+1, 2.4938e+1, 2.1528e+1, 2.0362e+1, |
---|
| 418 | 1.9474e+1, 1.8238e+1, 1.7898e+1, 1.9174e+1, 1.9023e+1, |
---|
| 419 | 1.8194e+1, 1.8504e+1, 1.8955e+1, 1.4276e+1, 1.7558e+1, |
---|
| 420 | 1.8651e+1, 1.7984e+1, 1.6793e+1, 1.5469e+1, 1.4143e+1, |
---|
| 421 | 1.3149e+1, 1.2255e+1, 9.2352e+0, 8.6067e+0, 9.7460e+0, |
---|
| 422 | 1.1749e+1, 1.3281e+1, 1.4326e+1, 1.4920e+1, 1.5157e+1, |
---|
| 423 | 1.5131e+1, 1.9489e+1, 2.3649e+1, 2.4686e+1, 2.4760e+1, |
---|
| 424 | 2.1519e+1, 2.0099e+1, 1.8746e+1, 1.5943e+1, 1.4880e+1}; |
---|
| 425 | |
---|
| 426 | G4double RA3[ntot] = {0.0e0, 2.2648e+0, 1.0579e+3, 8.6177e+2, 2.4422e+2, |
---|
| 427 | 7.8788e+1, 3.8293e+1, 1.2564e+1, 6.9091e+0, 3.7926e+0, |
---|
| 428 | 0.0000e+0, 0.0000e+0, 1.6759e-9, 1.3026e+1, 3.0569e+0, |
---|
| 429 | 1.5521e+2, 1.2815e+2, 4.7378e+1, 9.2802e+2, 4.7508e+2, |
---|
| 430 | 3.6612e+2, 2.7582e+2, 2.1008e+2, 1.5903e+2, 1.2322e+2, |
---|
| 431 | 9.2898e+1, 7.1345e+1, 5.1651e+1, 3.8474e+1, 2.7410e+1, |
---|
| 432 | 1.9126e+1, 1.0889e+1, 5.3479e+0, 8.2223e+0, 5.0837e+0, |
---|
| 433 | 2.8905e+0, 2.7457e+0, 6.7082e-1, 0.0000e+0, 0.0000e+0, |
---|
| 434 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 435 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 436 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 437 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 1.7264e-1, 2.7322e-1, |
---|
| 438 | 3.9444e-1, 4.5648e-1, 6.2286e-1, 7.2468e-1, 8.4296e-1, |
---|
| 439 | 1.1698e+0, 1.2994e+0, 1.4295e+0, 0.0000e+0, 8.1570e-1, |
---|
| 440 | 6.9349e-1, 4.9536e-1, 3.1211e-1, 1.5931e-1, 2.9512e-2, |
---|
| 441 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 442 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 443 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
| 444 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0}; |
---|
| 445 | |
---|
| 446 | G4double RA4[ntot] = {1.1055e1,6.3519e0,4.7367e+1, 3.9402e+1, 2.2896e+1, |
---|
| 447 | 1.3979e+1, 1.0766e+1, 6.5252e+0, 5.1631e+0, 4.0524e+0, |
---|
| 448 | 2.7145e+1, 1.8724e+1, 1.4782e+1, 1.1608e+1, 9.7750e+0, |
---|
| 449 | 1.6170e+1, 1.5249e+1, 9.1916e+0, 5.4499e+1, 4.1381e+1, |
---|
| 450 | 3.7395e+1, 3.3815e+1, 3.1135e+1, 2.8273e+1, 2.6140e+1, |
---|
| 451 | 2.3948e+1, 2.2406e+1, 2.0484e+1, 1.8453e+1, 1.7386e+1, |
---|
| 452 | 1.7301e+1, 1.5388e+1, 1.3411e+1, 1.2668e+1, 1.1133e+1, |
---|
| 453 | 9.8081e+0, 1.3031e+1, 1.4143e+1, 1.3815e+1, 1.2077e+1, |
---|
| 454 | 1.0033e+1, 9.2549e+0, 9.5338e+0, 7.9076e+0, 7.3263e+0, |
---|
| 455 | 5.9996e+0, 6.0087e+0, 6.2660e+0, 6.7914e+0, 7.1501e+0, |
---|
| 456 | 7.3367e+0, 7.3729e+0, 7.3508e+0, 7.1465e+0, 8.2731e+0, |
---|
| 457 | 9.3745e+0, 9.3508e+0, 8.9897e+0, 8.4566e+0, 8.2690e+0, |
---|
| 458 | 8.1398e+0, 7.9183e+0, 7.9123e+0, 8.1677e+0, 8.1871e+0, |
---|
| 459 | 8.1766e+0, 8.2881e+0, 8.4227e+0, 7.0273e+0, 8.0002e+0, |
---|
| 460 | 8.1440e+0, 7.9104e+0, 7.5685e+0, 7.1970e+0, 6.8184e+0, |
---|
| 461 | 6.5469e+0, 6.3056e+0, 5.4844e+0, 5.2832e+0, 5.5889e+0, |
---|
| 462 | 6.0919e+0, 6.4340e+0, 6.6426e+0, 6.7428e+0, 6.7636e+0, |
---|
| 463 | 6.7281e+0, 7.5729e+0, 8.2808e+0, 8.4400e+0, 8.4220e+0, |
---|
| 464 | 7.8662e+0, 7.5993e+0, 7.3353e+0, 6.7829e+0, 6.5520e+0}; |
---|
| 465 | |
---|
| 466 | G4double RA5[ntot] = {0.0e0, 4.9828e+0, 5.5674e+1, 3.0902e+1, 1.1496e+1, |
---|
| 467 | 4.8936e+0, 2.5506e+0, 1.2236e+0, 7.4698e-1, 4.7042e-1, |
---|
| 468 | 4.7809e+0, 4.6315e+0, 4.3677e+0, 4.9269e+0, 2.6033e+0, |
---|
| 469 | 9.6229e+0, 7.2592e+0, 4.1634e+0, 1.3999e+1, 8.6975e+0, |
---|
| 470 | 6.9630e+0, 5.4681e+0, 4.2653e+0, 3.2848e+0, 2.7354e+0, |
---|
| 471 | 2.1617e+0, 1.7030e+0, 1.2826e+0, 9.7080e-1, 7.2227e-1, |
---|
| 472 | 5.0874e-1, 3.1402e-1, 1.6360e-1, 3.2918e-1, 2.3570e-1, |
---|
| 473 | 1.5868e-1, 1.5146e-1, 9.7662e-2, 7.3151e-2, 6.4206e-2, |
---|
| 474 | 4.8945e-2, 4.3189e-2, 4.4368e-2, 3.3976e-2, 3.0466e-2, |
---|
| 475 | 2.4477e-2, 3.7202e-2, 3.7093e-2, 3.8161e-2, 3.8576e-2, |
---|
| 476 | 3.8403e-2, 3.7806e-2, 3.4958e-2, 3.6029e-2, 4.3087e-2, |
---|
| 477 | 4.7069e-2, 4.6452e-2, 4.2486e-2, 4.1517e-2, 4.1691e-2, |
---|
| 478 | 4.2813e-2, 4.2294e-2, 4.5287e-2, 4.8462e-2, 4.9726e-2, |
---|
| 479 | 5.5097e-2, 5.6568e-2, 5.8069e-2, 1.2270e-2, 3.8006e-2, |
---|
| 480 | 3.5048e-2, 3.0050e-2, 2.5069e-2, 2.0485e-2, 1.6151e-2, |
---|
| 481 | 1.4631e-2, 1.4034e-2, 1.1978e-2, 1.1522e-2, 1.2375e-2, |
---|
| 482 | 1.3805e-2, 1.4954e-2, 1.5832e-2, 1.6467e-2, 1.6896e-2, |
---|
| 483 | 1.7166e-2, 1.9954e-2, 2.2497e-2, 2.1942e-2, 2.1965e-2, |
---|
| 484 | 2.0005e-2, 1.8927e-2, 1.8167e-2, 1.6314e-2, 1.5522e-2}; |
---|
| 485 | |
---|
| 486 | G4double x=std::sqrt(y); |
---|
| 487 | G4double gradx1=0.0; |
---|
| 488 | G4double fa=0.0; |
---|
| 489 | |
---|
| 490 | G4int nElements = material->GetNumberOfElements(); |
---|
| 491 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
| 492 | const G4int* stechiometric = material->GetAtomsVector(); |
---|
| 493 | const G4double* vector_of_atoms = material->GetVecNbOfAtomsPerVolume(); |
---|
| 494 | const G4double tot_atoms = material->GetTotNbOfAtomsPerVolume(); |
---|
| 495 | for (G4int i=0;i<nElements;i++){ |
---|
| 496 | G4int Z = (G4int) (*elementVector)[i]->GetZ(); |
---|
| 497 | if (Z>ntot) Z=95; |
---|
| 498 | G4double denomin = 1+y*(RA4[Z-1]+y*RA5[Z-1]); |
---|
| 499 | fa=Z*(1+y*(RA1[Z-1]+x*(RA2[Z-1]+x*RA3[Z-1])))/(denomin*denomin); |
---|
| 500 | G4bool a = ((Z>10) && (fa>2.0)); |
---|
| 501 | if (a) |
---|
| 502 | { |
---|
| 503 | G4double Pa,Pg,Pq,fb; |
---|
| 504 | G4double k1=0.3125; |
---|
| 505 | G4double k2=2.426311e-02; |
---|
| 506 | Pa=(Z-k1)*fine_structure_const; |
---|
| 507 | Pg=std::sqrt(1-(Pa*Pa)); |
---|
| 508 | Pq=k2*x/Pa; |
---|
| 509 | fb=std::sin(2*Pg*std::atan(Pq))/(Pg*Pq*std::pow((1+Pq*Pq),Pg)); |
---|
| 510 | fa=std::max(fa,fb); |
---|
| 511 | } |
---|
| 512 | if (stechiometric && stechiometric[i]!=0) |
---|
| 513 | { |
---|
| 514 | gradx1=gradx1+stechiometric[i]*(fa*fa); //sum on the molecule |
---|
| 515 | } |
---|
| 516 | else |
---|
| 517 | { |
---|
| 518 | gradx1=gradx1+(vector_of_atoms[i]/tot_atoms)*(fa*fa); //weighted mean |
---|
| 519 | } |
---|
| 520 | } |
---|
| 521 | return gradx1; |
---|
| 522 | } |
---|
| 523 | |
---|
| 524 | |
---|
| 525 | G4double G4PenelopeRayleigh::DifferentialCrossSection(G4double x) |
---|
| 526 | { |
---|
| 527 | G4double x2=facte*(1-x); |
---|
| 528 | G4double gradx = (1+x*x)*MolecularFormFactor(x2); |
---|
| 529 | return gradx; |
---|
| 530 | } |
---|
| 531 | |
---|
| 532 | |
---|