00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050 #include "G4PenelopeCompton.hh"
00051 #include "Randomize.hh"
00052 #include "G4ParticleDefinition.hh"
00053 #include "G4Track.hh"
00054 #include "G4Step.hh"
00055 #include "G4ForceCondition.hh"
00056 #include "G4Gamma.hh"
00057 #include "G4Electron.hh"
00058 #include "G4DynamicParticle.hh"
00059 #include "G4VParticleChange.hh"
00060 #include "G4ThreeVector.hh"
00061 #include "G4EnergyLossTables.hh"
00062 #include "G4VCrossSectionHandler.hh"
00063 #include "G4CrossSectionHandler.hh"
00064 #include "G4VEMDataSet.hh"
00065 #include "G4EMDataSet.hh"
00066 #include "G4CompositeEMDataSet.hh"
00067 #include "G4VDataSetAlgorithm.hh"
00068 #include "G4LogLogInterpolation.hh"
00069 #include "G4VRangeTest.hh"
00070 #include "G4RangeTest.hh"
00071 #include "G4ProductionCutsTable.hh"
00072 #include "G4AtomicTransitionManager.hh"
00073 #include "G4AtomicShell.hh"
00074 #include "G4AtomicDeexcitation.hh"
00075 #include "G4PenelopeIntegrator.hh"
00076 #include "G4MaterialCutsCouple.hh"
00077
00078
00079 G4PenelopeCompton::G4PenelopeCompton(const G4String& processName)
00080 : G4VDiscreteProcess(processName),
00081 lowEnergyLimit(250*eV),
00082 highEnergyLimit(100*GeV),
00083 intrinsicLowEnergyLimit(10*eV),
00084 intrinsicHighEnergyLimit(100*GeV),
00085 energyForIntegration(0.0),
00086 ZForIntegration(1),
00087 nBins(200),
00088 cutForLowEnergySecondaryPhotons(250.0*eV)
00089 {
00090 if (lowEnergyLimit < intrinsicLowEnergyLimit ||
00091 highEnergyLimit > intrinsicHighEnergyLimit)
00092 {
00093 G4Exception("G4PenelopeCompton::G4PenelopeCompton - energy outside intrinsic process validity range");
00094 }
00095
00096 meanFreePathTable = 0;
00097 ionizationEnergy = new std::map<G4int,G4DataVector*>;
00098 hartreeFunction = new std::map<G4int,G4DataVector*>;
00099 occupationNumber = new std::map<G4int,G4DataVector*>;
00100
00101 rangeTest = new G4RangeTest;
00102
00103 ReadData();
00104
00105 if (verboseLevel > 0)
00106 {
00107 G4cout << GetProcessName() << " is created " << G4endl
00108 << "Energy range: "
00109 << lowEnergyLimit / keV << " keV - "
00110 << highEnergyLimit / GeV << " GeV"
00111 << G4endl;
00112 }
00113 }
00114
00115 G4PenelopeCompton::~G4PenelopeCompton()
00116 {
00117 delete meanFreePathTable;
00118 delete rangeTest;
00119
00120 for (size_t i1=0;i1<matCrossSections->size();i1++)
00121 {
00122 delete (*matCrossSections)[i1];
00123 }
00124
00125 delete matCrossSections;
00126
00127 for (G4int Z=1;Z<100;Z++)
00128 {
00129 if (ionizationEnergy->count(Z)) delete (ionizationEnergy->find(Z)->second);
00130 if (hartreeFunction->count(Z)) delete (hartreeFunction->find(Z)->second);
00131 if (occupationNumber->count(Z)) delete (occupationNumber->find(Z)->second);
00132 }
00133 delete ionizationEnergy;
00134 delete hartreeFunction;
00135 delete occupationNumber;
00136 }
00137
00138 void G4PenelopeCompton::BuildPhysicsTable(const G4ParticleDefinition& )
00139 {
00140 G4DataVector energyVector;
00141 G4double dBin = std::log10(highEnergyLimit/lowEnergyLimit)/nBins;
00142 G4int i;
00143 for (i=0;i<nBins;i++)
00144 {
00145 energyVector.push_back(std::pow(10.,std::log10(lowEnergyLimit)+i*dBin));
00146 }
00147
00148 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
00149 G4int nMaterials = G4Material::GetNumberOfMaterials();
00150 G4VDataSetAlgorithm* algo = new G4LogLogInterpolation();
00151
00152 size_t nOfBins = energyVector.size();
00153 size_t bin=0;
00154
00155 G4DataVector* energies;
00156 G4DataVector* data;
00157
00158 matCrossSections = new std::vector<G4VEMDataSet*>;
00159
00160
00161 G4int m;
00162 for (m=0; m<nMaterials; m++)
00163 {
00164 const G4Material* material= (*materialTable)[m];
00165 G4int nElements = material->GetNumberOfElements();
00166 const G4ElementVector* elementVector = material->GetElementVector();
00167 const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
00168
00169 G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
00170
00171 for (i=0; i<nElements; i++) {
00172
00173 G4int Z = (G4int) (*elementVector)[i]->GetZ();
00174 G4double density = nAtomsPerVolume[i];
00175 G4double cross=0.0;
00176 energies = new G4DataVector;
00177 data = new G4DataVector;
00178
00179
00180 for (bin=0; bin<nOfBins; bin++)
00181 {
00182 G4double e = energyVector[bin];
00183 energies->push_back(e);
00184 cross = density * CrossSection(e,Z);
00185 data->push_back(cross);
00186 }
00187
00188 G4VEMDataSet* elSet = new G4EMDataSet(i,energies,data,algo,1.,1.);
00189 setForMat->AddComponent(elSet);
00190 }
00191
00192 matCrossSections->push_back(setForMat);
00193 }
00194
00195
00196
00197 G4double matCS = 0.0;
00198 G4VEMDataSet* matCrossSet = new G4CompositeEMDataSet(algo,1.,1.);
00199 G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo,1.,1.);
00200
00201
00202 for (m=0; m<nMaterials; m++)
00203 {
00204 energies = new G4DataVector;
00205 data = new G4DataVector;
00206 const G4Material* material= (*materialTable)[m];
00207 material= (*materialTable)[m];
00208 for (bin=0; bin<nOfBins; bin++)
00209 {
00210 G4double energy = energyVector[bin];
00211 energies->push_back(energy);
00212 matCrossSet = (*matCrossSections)[m];
00213 matCS = 0.0;
00214 G4int nElm = matCrossSet->NumberOfComponents();
00215 for(G4int j=0; j<nElm; j++) {
00216 matCS += matCrossSet->GetComponent(j)->FindValue(energy);
00217 }
00218 if (matCS > 0.)
00219 {
00220 data->push_back(1./matCS);
00221 }
00222 else
00223 {
00224 data->push_back(DBL_MAX);
00225 }
00226 }
00227 G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.);
00228 materialSet->AddComponent(dataSet);
00229 }
00230 meanFreePathTable = materialSet;
00231 }
00232
00233 G4VParticleChange* G4PenelopeCompton::PostStepDoIt(const G4Track& aTrack,
00234 const G4Step& aStep)
00235 {
00236
00237
00238 aParticleChange.Initialize(aTrack);
00239
00240
00241 const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
00242 G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
00243
00244 if (photonEnergy0 <= lowEnergyLimit)
00245 {
00246 aParticleChange.ProposeTrackStatus(fStopAndKill);
00247 aParticleChange.ProposeEnergy(0.);
00248 aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
00249 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
00250 }
00251
00252 G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
00253
00254 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
00255 const G4Material* material = couple->GetMaterial();
00256 G4int Z = SelectRandomAtomForCompton(material,photonEnergy0);
00257 const G4int nmax = 64;
00258 G4double rn[nmax],pac[nmax];
00259
00260 G4double ki,ki1,ki2,ki3,taumin,a1,a2;
00261 G4double tau,TST;
00262 G4double S=0.0;
00263 G4double epsilon,cosTheta;
00264 G4double harFunc = 0.0;
00265 G4int occupNb= 0;
00266 G4double ionEnergy=0.0;
00267 G4int nosc = occupationNumber->find(Z)->second->size();
00268 G4int iosc = nosc;
00269 ki = photonEnergy0/electron_mass_c2;
00270 ki2 = 2*ki+1.0;
00271 ki3 = ki*ki;
00272 ki1 = ki3-ki2-1.0;
00273 taumin = 1.0/ki2;
00274 a1 = std::log(ki2);
00275 a2 = a1+2.0*ki*(1.0+ki)/(ki2*ki2);
00276 if (photonEnergy0 > 5*MeV)
00277 {
00278 do{
00279 do{
00280 if ((a2*G4UniformRand()) < a1)
00281 {
00282 tau = std::pow(taumin,G4UniformRand());
00283 }
00284 else
00285 {
00286 tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
00287 }
00288
00289 TST = (1+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau));
00290 }while (G4UniformRand()> TST);
00291 epsilon=tau;
00292 cosTheta = 1.0 - (1.0-tau)/(ki*tau);
00293
00294 TST = Z*G4UniformRand();
00295 iosc = nosc;
00296 S=0.0;
00297 for (G4int j=0;j<nosc;j++)
00298 {
00299 occupNb = (G4int) (*(occupationNumber->find(Z)->second))[j];
00300 S = S + occupNb;
00301 if (S > TST) iosc = j;
00302 if (S > TST) break;
00303 }
00304 ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc];
00305 }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0);
00306 }
00307
00308 else
00309 {
00310
00311 G4double s0=0.0;
00312 G4double pzomc=0.0,rni=0.0;
00313 G4double aux=0.0;
00314 for (G4int i=0;i<nosc;i++){
00315 ionEnergy = (*(ionizationEnergy->find(Z)->second))[i];
00316 if (photonEnergy0 > ionEnergy)
00317 {
00318 G4double aux = photonEnergy0*(photonEnergy0-ionEnergy)*2.0;
00319 harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const;
00320 occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
00321 pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/
00322 (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy));
00323 if (pzomc > 0)
00324 {
00325 rni = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
00326 }
00327 else
00328 {
00329 rni = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
00330 }
00331 s0 = s0 + occupNb*rni;
00332 }
00333 }
00334
00335
00336 G4double cdt1;
00337 do
00338 {
00339 if ((G4UniformRand()*a2) < a1)
00340 {
00341 tau = std::pow(taumin,G4UniformRand());
00342 }
00343 else
00344 {
00345 tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
00346 }
00347 cdt1 = (1.0-tau)/(ki*tau);
00348 S=0.0;
00349
00350 for (G4int i=0;i<nosc;i++){
00351 ionEnergy = (*(ionizationEnergy->find(Z)->second))[i];
00352 if (photonEnergy0 > ionEnergy)
00353 {
00354 aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1;
00355 harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const;
00356 occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
00357 pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/
00358 (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy));
00359 if (pzomc > 0)
00360 {
00361 rn[i] = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
00362 }
00363 else
00364 {
00365 rn[i] = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
00366 }
00367 S = S + occupNb*rn[i];
00368 pac[i] = S;
00369 }
00370 else
00371 {
00372 pac[i] = S-(1e-06);
00373 }
00374 }
00375
00376 TST = S*(1.0+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau));
00377 }while ((G4UniformRand()*s0) > TST);
00378
00379
00380 cosTheta = 1.0 - cdt1;
00381 G4double fpzmax=0.0,fpz=0.0;
00382 G4double A=0.0;
00383 do
00384 {
00385 do
00386 {
00387 TST =S*G4UniformRand();
00388 iosc=nosc;
00389 for (G4int i=0;i<nosc;i++){
00390 if (pac[i]>TST) iosc = i;
00391 if (pac[i]>TST) break;
00392 }
00393 A = G4UniformRand()*rn[iosc];
00394 harFunc = (*(hartreeFunction->find(Z)->second))[iosc]/fine_structure_const;
00395 occupNb = (G4int) (*(occupationNumber->find(Z)->second))[iosc];
00396 if (A < 0.5) {
00397 pzomc = (std::sqrt(0.5)-std::sqrt(0.5-std::log(2.0*A)))/
00398 (std::sqrt(2.0)*harFunc);
00399 }
00400 else
00401 {
00402 pzomc = (std::sqrt(0.5-std::log(2.0-2.0*A))-std::sqrt(0.5))/
00403 (std::sqrt(2.0)*harFunc);
00404 }
00405 } while (pzomc < -1);
00406
00407 G4double XQC = 1.0+tau*(tau-2.0*cosTheta);
00408 G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC);
00409 if (AF > 0) {
00410 fpzmax = 1.0+AF*0.2;
00411 }
00412 else
00413 {
00414 fpzmax = 1.0-AF*0.2;
00415 }
00416 fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2);
00417 }while ((fpzmax*G4UniformRand())>fpz);
00418
00419
00420 G4double T = pzomc*pzomc;
00421 G4double b1 = 1.0-T*tau*tau;
00422 G4double b2 = 1.0-T*tau*cosTheta;
00423 if (pzomc > 0.0)
00424 {
00425 epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
00426 }
00427 else
00428 {
00429 epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
00430 }
00431 }
00432
00433
00434 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
00435 G4double phi = twopi * G4UniformRand() ;
00436 G4double dirx = sinTheta * std::cos(phi);
00437 G4double diry = sinTheta * std::sin(phi);
00438 G4double dirz = cosTheta ;
00439
00440
00441
00442 G4ThreeVector photonDirection1(dirx,diry,dirz);
00443 photonDirection1.rotateUz(photonDirection0);
00444 aParticleChange.ProposeMomentumDirection(photonDirection1) ;
00445 G4double photonEnergy1 = epsilon * photonEnergy0;
00446
00447 if (photonEnergy1 > 0.)
00448 {
00449 aParticleChange.ProposeEnergy(photonEnergy1) ;
00450 }
00451 else
00452 {
00453 aParticleChange.ProposeEnergy(0.) ;
00454 aParticleChange.ProposeTrackStatus(fStopAndKill);
00455 }
00456
00457
00458
00459
00460
00461 G4double diffEnergy = photonEnergy0*(1-epsilon);
00462 ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc];
00463
00464 G4double Q2 = photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta);
00465 G4double cosThetaE;
00466 if (Q2 > 1.0e-12)
00467 {
00468 cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2);
00469 }
00470 else
00471 {
00472 cosThetaE = 1.0;
00473 }
00474 G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE);
00475
00476
00477
00478 const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
00479 const G4AtomicShell* shell = transitionManager->Shell(Z,iosc);
00480 G4double bindingEnergy = shell->BindingEnergy();
00481 G4int shellId = shell->ShellId();
00482
00483 ionEnergy = std::max(bindingEnergy,ionEnergy);
00484 G4double eKineticEnergy = diffEnergy - ionEnergy;
00485
00486 size_t nTotPhotons=0;
00487 G4int nPhotons=0;
00488
00489 const G4ProductionCutsTable* theCoupleTable=
00490 G4ProductionCutsTable::GetProductionCutsTable();
00491 size_t indx = couple->GetIndex();
00492 G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx];
00493 cutg = std::min(cutForLowEnergySecondaryPhotons,cutg);
00494
00495 G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[indx];
00496 cute = std::min(cutForLowEnergySecondaryPhotons,cute);
00497
00498 std::vector<G4DynamicParticle*>* photonVector=0;
00499 G4DynamicParticle* aPhoton;
00500 G4AtomicDeexcitation deexcitationManager;
00501
00502 if (Z>5 && (ionEnergy > cutg || ionEnergy > cute))
00503 {
00504 photonVector = deexcitationManager.GenerateParticles(Z,shellId);
00505 nTotPhotons = photonVector->size();
00506 for (size_t k=0;k<nTotPhotons;k++){
00507 aPhoton = (*photonVector)[k];
00508 if (aPhoton)
00509 {
00510 G4double itsCut = cutg;
00511 if (aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cute;
00512 G4double itsEnergy = aPhoton->GetKineticEnergy();
00513 if (itsEnergy > itsCut && itsEnergy <= ionEnergy)
00514 {
00515 nPhotons++;
00516 ionEnergy -= itsEnergy;
00517 }
00518 else
00519 {
00520 delete aPhoton;
00521 (*photonVector)[k]=0;
00522 }
00523 }
00524 }
00525 }
00526 G4double energyDeposit =ionEnergy;
00527 G4int nbOfSecondaries=nPhotons;
00528
00529
00530 G4double safety = aStep.GetPostStepPoint()->GetSafety();
00531 G4DynamicParticle* electron = 0;
00532 if (rangeTest->Escape(G4Electron::Electron(),couple,eKineticEnergy,safety))
00533 {
00534 G4double xEl = sinThetaE * std::cos(phi+pi);
00535 G4double yEl = sinThetaE * std::sin(phi+pi);
00536 G4double zEl = cosThetaE;
00537 G4ThreeVector eDirection(xEl,yEl,zEl);
00538 eDirection.rotateUz(photonDirection0);
00539 electron = new G4DynamicParticle (G4Electron::Electron(),
00540 eDirection,eKineticEnergy) ;
00541 nbOfSecondaries++;
00542 }
00543 else
00544 {
00545
00546 energyDeposit += eKineticEnergy;
00547 }
00548
00549 aParticleChange.SetNumberOfSecondaries(nbOfSecondaries);
00550 if (electron) aParticleChange.AddSecondary(electron);
00551 for (size_t ll=0;ll<nTotPhotons;ll++)
00552 {
00553 aPhoton = (*photonVector)[ll];
00554 if (aPhoton) aParticleChange.AddSecondary(aPhoton);
00555 }
00556 delete photonVector;
00557 if (energyDeposit < 0)
00558 {
00559 G4cout << "WARNING-"
00560 << "G4PenelopeCompton::PostStepDoIt - Negative energy deposit"
00561 << G4endl;
00562 energyDeposit=0;
00563 }
00564 aParticleChange.ProposeLocalEnergyDeposit(energyDeposit);
00565
00566
00567 return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
00568 }
00569
00570 G4bool G4PenelopeCompton::IsApplicable(const G4ParticleDefinition& particle)
00571 {
00572 return ( &particle == G4Gamma::Gamma() );
00573 }
00574
00575 G4double G4PenelopeCompton::GetMeanFreePath(const G4Track& track,
00576 G4double,
00577 G4ForceCondition*)
00578 {
00579 const G4DynamicParticle* photon = track.GetDynamicParticle();
00580 G4double energy = photon->GetKineticEnergy();
00581 G4Material* material = track.GetMaterial();
00582 size_t materialIndex = material->GetIndex();
00583
00584 G4double meanFreePath;
00585 if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
00586 else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
00587 else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
00588 return meanFreePath;
00589 }
00590
00591
00592 void G4PenelopeCompton::ReadData()
00593 {
00594 char* path = getenv("G4LEDATA");
00595 if (!path)
00596 {
00597 G4String excep = "G4PenelopeCompton - G4LEDATA environment variable not set!";
00598 G4Exception(excep);
00599 }
00600 G4String pathString(path);
00601 G4String pathFile = pathString + "/penelope/compton-pen.dat";
00602 std::ifstream file(pathFile);
00603 std::filebuf* lsdp = file.rdbuf();
00604
00605 if (!(lsdp->is_open()))
00606 {
00607 G4String excep = "G4PenelopeCompton - data file " + pathFile + " not found!";
00608 G4Exception(excep);
00609 }
00610
00611 G4int k1,test,test1;
00612 G4double a1,a2;
00613 G4int Z=1,nLevels=0;
00614 G4DataVector* f;
00615 G4DataVector* u;
00616 G4DataVector* j;
00617
00618 do{
00619 f = new G4DataVector;
00620 u = new G4DataVector;
00621 j = new G4DataVector;
00622 file >> Z >> nLevels;
00623 for (G4int h=0;h<nLevels;h++){
00624 file >> k1 >> a1 >> a2;
00625 f->push_back((G4double) k1);
00626 u->push_back(a1);
00627 j->push_back(a2);
00628 }
00629 ionizationEnergy->insert(std::make_pair(Z,u));
00630 hartreeFunction->insert(std::make_pair(Z,j));
00631 occupationNumber->insert(std::make_pair(Z,f));
00632 file >> test >> test1;
00633 if (test > 0) {
00634 G4String excep = "G4PenelopeCompton - data file corrupted!";
00635 G4Exception(excep);
00636 }
00637 }while (test != -2);
00638 }
00639
00640 G4double G4PenelopeCompton::CrossSection(G4double energy,G4int Z)
00641 {
00642 G4double cs=0.0;
00643 energyForIntegration=energy;
00644 ZForIntegration = Z;
00645 if (energy< 5*MeV)
00646 {
00647 G4PenelopeIntegrator<G4PenelopeCompton,G4double (G4PenelopeCompton::*)(G4double)> theIntegrator;
00648 cs = theIntegrator.Calculate(this,&G4PenelopeCompton::DifferentialCrossSection,-1.0,1.0,1e-05);
00649 }
00650 else
00651 {
00652 G4double ki=energy/electron_mass_c2;
00653 G4double ki3=ki*ki;
00654 G4double ki2=1.0+2*ki;
00655 G4double ki1=ki3-ki2-1.0;
00656 G4double t0=1.0/(ki2);
00657 G4double csl = 0.5*ki3*t0*t0+ki2*t0+ki1*std::log(t0)-(1.0/t0);
00658 G4int nosc = occupationNumber->find(Z)->second->size();
00659 for (G4int i=0;i<nosc;i++)
00660 {
00661 G4double ionEnergy = (*(ionizationEnergy->find(Z)->second))[i];
00662 G4double tau=(energy-ionEnergy)/energy;
00663 if (tau > t0)
00664 {
00665 G4double csu = 0.5*ki3*tau*tau+ki2*tau+ki1*std::log(tau)-(1.0/tau);
00666 G4int f = (G4int) (*(occupationNumber->find(Z)->second))[i];
00667 cs = cs + f*(csu-csl);
00668 }
00669 }
00670 cs=pi*classic_electr_radius*classic_electr_radius*cs/(ki*ki3);
00671 }
00672 return cs;
00673 }
00674
00675
00676 G4double G4PenelopeCompton::DifferentialCrossSection(G4double cosTheta)
00677 {
00678 const G4double k2 = std::sqrt(2.0);
00679 const G4double k1 = std::sqrt(0.5);
00680 const G4double k12 = 0.5;
00681 G4double cdt1 = 1.0-cosTheta;
00682 G4double energy = energyForIntegration;
00683 G4int Z = ZForIntegration;
00684 G4double ionEnergy=0.0,Pzimax=0.0,XKN=0.0;
00685 G4double diffCS=0.0;
00686 G4double x=0.0,siap=0.0;
00687 G4double harFunc=0.0;
00688 G4int occupNb;
00689
00690 G4double EOEC = 1.0+(energy/electron_mass_c2)*cdt1;
00691 G4double ECOE = 1.0/EOEC;
00692
00693 G4double sia = 0.0;
00694 G4int nosc = occupationNumber->find(Z)->second->size();
00695 for (G4int i=0;i<nosc;i++){
00696 ionEnergy = (*(ionizationEnergy->find(Z)->second))[i];
00697
00698 if (energy > ionEnergy)
00699 {
00700 G4double aux = energy * (energy-ionEnergy)*cdt1;
00701 Pzimax = (aux - electron_mass_c2*ionEnergy)/(electron_mass_c2*std::sqrt(2*aux+ionEnergy*ionEnergy));
00702 harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const;
00703 occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
00704 x = harFunc*Pzimax;
00705 if (x > 0)
00706 {
00707 siap = 1.0-0.5*std::exp(k12-(k1+k2*x)*(k1+k2*x));
00708 }
00709 else
00710 {
00711 siap = 0.5*std::exp(k12-(k1-k2*x)*(k1-k2*x));
00712 }
00713 sia = sia + occupNb*siap;
00714 }
00715 }
00716 XKN = EOEC+ECOE-1+cosTheta*cosTheta;
00717 diffCS = pi*classic_electr_radius*classic_electr_radius*ECOE*ECOE*XKN*sia;
00718 return diffCS;
00719 }
00720
00721 G4int G4PenelopeCompton::SelectRandomAtomForCompton(const G4Material* material,G4double energy) const
00722 {
00723 G4int nElements = material->GetNumberOfElements();
00724
00725 if (nElements == 1)
00726 {
00727 G4int Z = (G4int) material->GetZ();
00728 return Z;
00729 }
00730
00731
00732 const G4ElementVector* elementVector = material->GetElementVector();
00733 size_t materialIndex = material->GetIndex();
00734
00735 G4VEMDataSet* materialSet = (*matCrossSections)[materialIndex];
00736 G4double materialCrossSection0 = 0.0;
00737 G4DataVector cross;
00738 cross.clear();
00739 G4int i;
00740 for (i=0;i<nElements;i++)
00741 {
00742 G4double cr = (materialSet->GetComponent(i))->FindValue(energy);
00743 materialCrossSection0 += cr;
00744 cross.push_back(materialCrossSection0);
00745 }
00746
00747 G4double random = G4UniformRand()*materialCrossSection0;
00748 for (i=0;i<nElements;i++)
00749 {
00750 if (random <= cross[i]) return (G4int) (*elementVector)[i]->GetZ();
00751 }
00752
00753 return 0;
00754 }
00755