62 :
G4VEmModel(nam),fParticleChange(0),isInitialised(false),fAtomDeexcitation(0),
65 fIntrinsicLowEnergyLimit = 100.0*eV;
66 fIntrinsicHighEnergyLimit = 100.0*GeV;
88 std::map <G4int,G4PhysicsTable*>::iterator i;
91 for (i=logAtomicShellXS->begin();i != logAtomicShellXS->end();i++)
98 delete logAtomicShellXS;
106 if (verboseLevel > 3)
107 G4cout <<
"Calling G4PenelopePhotoElectricModel::Initialise()" <<
G4endl;
110 if (!logAtomicShellXS)
111 logAtomicShellXS =
new std::map<G4int,G4PhysicsTable*>;
117 if (!fAtomDeexcitation)
120 G4cout <<
"WARNING from G4PenelopePhotoElectricModel " <<
G4endl;
121 G4cout <<
"Atomic de-excitation module is not instantiated, so there will not be ";
123 G4cout <<
"Please make sure this is intended" <<
G4endl;
126 if (verboseLevel > 0) {
127 G4cout <<
"Penelope Photo-Electric model v2008 is initialized " <<
G4endl
133 if(isInitialised)
return;
135 isInitialised =
true;
151 if (verboseLevel > 3)
152 G4cout <<
"Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" <<
G4endl;
157 if (!logAtomicShellXS->count(iZ))
160 if (!logAtomicShellXS->count(iZ))
161 G4Exception(
"G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
163 "Unable to retrieve the shell cross section table");
172 G4Exception(
"G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
174 "Unable to retrieve the total cross section table");
179 cross = std::exp(logXS);
181 if (verboseLevel > 2)
182 G4cout <<
"Photoelectric cross section at " << energy/MeV <<
" MeV for Z=" << Z <<
183 " = " << cross/barn <<
" barn" <<
G4endl;
209 if (verboseLevel > 3)
210 G4cout <<
"Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" <<
G4endl;
218 if (photonEnergy <= fIntrinsicLowEnergyLimit)
227 if (verboseLevel > 2)
234 if (verboseLevel > 2)
243 size_t shellIndex = SelectRandomShell(Z,photonEnergy);
245 if (verboseLevel > 2)
246 G4cout <<
"Selected shell " << shellIndex <<
" of element " << anElement->
GetName() <<
G4endl;
257 size_t numberOfShells = (size_t) transitionManager->
NumberOfShells(Z);
258 if (shellIndex >= numberOfShells)
259 shellIndex = numberOfShells-1;
271 bindingEnergy = 0.*eV;
278 G4double eKineticEnergy = photonEnergy - bindingEnergy;
282 if (eKineticEnergy > 0.)
286 cosTheta = SampleElectronDirection(eKineticEnergy);
287 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
289 G4double dirx = sinTheta * std::cos(phi);
290 G4double diry = sinTheta * std::sin(phi);
293 electronDirection.
rotateUz(photonDirection);
297 fvect->push_back(electron);
300 bindingEnergy = photonEnergy;
309 if (fAtomDeexcitation && shellIndex<9)
314 size_t nBefore = fvect->size();
316 size_t nAfter = fvect->size();
318 if (nAfter > nBefore)
320 for (
size_t j=nBefore;j<nAfter;j++)
322 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
323 bindingEnergy -= itsEnergy;
325 energyInFluorescence += itsEnergy;
327 energyInAuger += itsEnergy;
336 localEnergyDeposit += bindingEnergy;
338 if (localEnergyDeposit < 0)
341 <<
"G4PenelopePhotoElectricModel::SampleSecondaries() - Negative energy deposit"
343 localEnergyDeposit = 0;
348 if (verboseLevel > 1)
350 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
351 G4cout <<
"Energy balance from G4PenelopePhotoElectric" <<
G4endl;
352 G4cout <<
"Selected shell: " << WriteTargetShell(shellIndex) <<
" of element " <<
354 G4cout <<
"Incoming photon energy: " << photonEnergy/keV <<
" keV" <<
G4endl;
355 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
357 G4cout <<
"Outgoing electron " << eKineticEnergy/keV <<
" keV" <<
G4endl;
358 if (energyInFluorescence)
359 G4cout <<
"Fluorescence x-rays: " << energyInFluorescence/keV <<
" keV" <<
G4endl;
361 G4cout <<
"Auger electrons: " << energyInAuger/keV <<
" keV" <<
G4endl;
362 G4cout <<
"Local energy deposit " << localEnergyDeposit/keV <<
" keV" <<
G4endl;
363 G4cout <<
"Total final state: " <<
364 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
366 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
368 if (verboseLevel > 0)
371 std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger-photonEnergy);
372 if (energyDiff > 0.05*keV)
374 G4cout <<
"Warning from G4PenelopePhotoElectric: problem with energy conservation: " <<
375 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV
376 <<
" keV (final) vs. " <<
377 photonEnergy/keV <<
" keV (initial)" <<
G4endl;
378 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
379 G4cout <<
"Energy balance from G4PenelopePhotoElectric" <<
G4endl;
380 G4cout <<
"Selected shell: " << WriteTargetShell(shellIndex) <<
" of element " <<
382 G4cout <<
"Incoming photon energy: " << photonEnergy/keV <<
" keV" <<
G4endl;
383 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
385 G4cout <<
"Outgoing electron " << eKineticEnergy/keV <<
" keV" <<
G4endl;
386 if (energyInFluorescence)
387 G4cout <<
"Fluorescence x-rays: " << energyInFluorescence/keV <<
" keV" <<
G4endl;
389 G4cout <<
"Auger electrons: " << energyInAuger/keV <<
" keV" <<
G4endl;
390 G4cout <<
"Local energy deposit " << localEnergyDeposit/keV <<
" keV" <<
G4endl;
391 G4cout <<
"Total final state: " <<
392 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
394 G4cout <<
"-----------------------------------------------------------" <<
G4endl;
401G4double G4PenelopePhotoElectricModel::SampleElectronDirection(
G4double energy)
404 if (energy>1*GeV)
return costheta;
409 G4double gamma = 1.0 + energy/electron_mass_c2;
411 G4double beta = std::sqrt((gamma2-1.0)/gamma2);
416 G4double a1 = 0.5*beta*gamma*(gamma-1.0)*(gamma-2.0);
428 tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand);
429 gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam));
439void G4PenelopePhotoElectricModel::ReadDataFile(
G4int Z)
441 if (verboseLevel > 2)
443 G4cout <<
"G4PenelopePhotoElectricModel::ReadDataFile()" <<
G4endl;
444 G4cout <<
"Going to read PhotoElectric data files for Z=" << Z <<
G4endl;
447 char* path = getenv(
"G4LEDATA");
450 G4String excep =
"G4PenelopePhotoElectricModel - G4LEDATA environment variable not set!";
451 G4Exception(
"G4PenelopePhotoElectricModel::ReadDataFile()",
459 std::ostringstream ost;
461 ost << path <<
"/penelope/photoelectric/pdgph" << Z <<
".p08";
463 ost << path <<
"/penelope/photoelectric/pdgph0" << Z <<
".p08";
464 std::ifstream file(ost.str().c_str());
467 G4String excep =
"G4PenelopePhotoElectricModel - data file " +
G4String(ost.str()) +
" not found!";
468 G4Exception(
"G4PenelopePhotoElectricModel::ReadDataFile()",
475 while( getline(file, line) )
482 file.open(ost.str().c_str());
486 file >> readZ >> nShells;
488 if (verboseLevel > 3)
489 G4cout <<
"Element Z=" << Z <<
" , nShells = " << nShells <<
G4endl;
492 if (readZ != Z || nShells <= 0 || nShells > 50)
495 ed <<
"Corrupted data file for Z=" << Z <<
G4endl;
496 G4Exception(
"G4PenelopePhotoElectricModel::ReadDataFile()",
508 for (
size_t i=0;i<nShells+1;i++)
512 for (k=0;k<ndata && !file.eof();k++)
520 for (
size_t i=0;i<nShells+1;i++)
525 if (aValue < 1e-40*cm2)
527 theVec->
PutValue(k,logene,std::log(aValue));
531 if (verboseLevel > 2)
533 G4cout <<
"G4PenelopePhotoElectricModel: read " << k <<
" points for element Z = "
537 logAtomicShellXS->insert(std::make_pair(Z,thePhysicsTable));
545size_t G4PenelopePhotoElectricModel::SelectRandomShell(
G4int Z,
G4double energy)
547 G4double logEnergy = std::log(energy);
550 if (!logAtomicShellXS->count(Z))
553 ed <<
"Cannot find shell cross section data for Z=" << Z <<
G4endl;
554 G4Exception(
"G4PenelopePhotoElectricModel::SelectRandomShell()",
558 size_t shellIndex = 0;
567 tempVector->push_back(sum);
580 for (
size_t k=1;k<theTable->
entries();k++)
584 G4double partialXS = std::exp(logXSLocal);
586 tempVector->push_back(sum);
589 tempVector->push_back(totalXS);
601 size_t lowerBound = 0;
602 size_t upperBound = tempVector->size()-1;
603 while (lowerBound <= upperBound)
605 size_t midBin = (lowerBound + upperBound)/2;
606 if( random < (*tempVector)[midBin])
607 upperBound = midBin-1;
609 lowerBound = midBin+1;
612 shellIndex = upperBound;
623 if (!logAtomicShellXS->count(Z))
626 if (!logAtomicShellXS->count(Z))
629 ed <<
"Cannot find shell cross section data for Z=" << Z <<
G4endl;
630 G4Exception(
"G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
634 size_t nEntries = logAtomicShellXS->find(Z)->second->entries();
645 if (shellID >= entries)
647 G4cout <<
"Element Z=" << Z <<
" has data for " << entries <<
" shells only" <<
G4endl;
648 G4cout <<
"so shellID should be from 0 to " << entries-1 <<
G4endl;
658 G4Exception(
"G4PenelopePhotoElectricModel::GetShellCrossSection()",
660 "Unable to retrieve the total cross section table");
666 if (cross < 2e-40*cm2) cross = 0;
672G4String G4PenelopePhotoElectricModel::WriteTargetShell(
size_t shellID)
677 else if (shellID == 1)
679 else if (shellID == 2)
681 else if (shellID == 3)
683 else if (shellID == 4)
685 else if (shellID == 5)
687 else if (shellID == 6)
689 else if (shellID == 7)
691 else if (shellID == 8)
G4DLLIMPORT std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
G4double BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
G4int NumberOfShells(G4int Z) const
static G4AtomicTransitionManager * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Definition()
static G4Electron * Electron()
const G4String & GetName() const
static G4Gamma * GammaDefinition()
static G4Gamma * Definition()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4PenelopePhotoElectricModel(const G4ParticleDefinition *p=0, const G4String &processName="PenPhotoElec")
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4double GetShellCrossSection(G4int Z, size_t shellID, G4double energy)
size_t GetNumberOfShellXS(G4int)
G4ParticleChangeForGamma * fParticleChange
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual ~G4PenelopePhotoElectricModel()
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
void push_back(G4PhysicsVector *)
G4double Value(G4double theEnergy)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetDeexcitationFlag(G4bool val)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
std::ostringstream G4ExceptionDescription