52 fParticleChange =
nullptr;
54 lowEnergyLimit = 2.0*electron_mass_c2;
64 G4cout <<
"G4LivermoreNuclearGammaConversionModel is constructed " <<
G4endl;
73 for(
G4int i=0; i<maxZ; ++i) {
90 G4cout <<
"Calling Initialise() of G4LivermoreNuclearGammaConversionModel."
112 for(
G4int i=0; i<numOfCouples; ++i)
119 for (std::size_t j=0; j<nelm; ++j)
123 else if(
Z > maxZ) {
Z = maxZ; }
124 if(!data[
Z]) { ReadData(
Z, path); }
128 if(isInitialised) {
return; }
130 isInitialised =
true;
148 return lowEnergyLimit;
153void G4LivermoreNuclearGammaConversionModel::ReadData(
size_t Z,
const char* path)
155 if (verboseLevel > 1)
157 G4cout <<
"Calling ReadData() of G4LivermoreNuclearGammaConversionModel"
162 if(data[
Z]) {
return; }
164 const char* datadir = path;
171 G4Exception(
"G4LivermoreNuclearGammaConversionModel::ReadData()",
173 "Environment variable G4LEDATA not defined");
180 std::ostringstream ost;
181 ost << datadir <<
"/livermore/pairdata/pp-pair-cs-" <<
Z <<
".dat";
182 std::ifstream fin(ost.str().c_str());
187 ed <<
"G4LivermoreNuclearGammaConversionModel data file <" << ost.str().c_str()
188 <<
"> is not opened!" <<
G4endl;
189 G4Exception(
"G4LivermoreNuclearGammaConversionModel::ReadData()",
191 ed,
"G4LEDATA version should be G4EMLOW8.0 or later.");
197 if(verboseLevel > 3) {
G4cout <<
"File " << ost.str()
198 <<
" is opened by G4LivermoreNuclearGammaConversionModel" <<
G4endl;}
215 if (verboseLevel > 1)
217 G4cout <<
"Calling ComputeCrossSectionPerAtom() of G4LivermoreNuclearGammaConversionModel"
221 if (GammaEnergy < lowEnergyLimit) {
return 0.0; }
227 if(intZ < 1 || intZ > maxZ) {
return xs; }
237 if(!pv) {
return xs; }
240 xs = pv->
Value(GammaEnergy);
245 G4cout <<
"****** DEBUG: tcs value for Z=" <<
Z <<
" at energy (MeV)="
246 << GammaEnergy/MeV <<
G4endl;
247 G4cout <<
" cs (Geant4 internal unit)=" << xs <<
G4endl;
248 G4cout <<
" -> first cs value in EADL data file (iu) =" << (*pv)[0] <<
G4endl;
249 G4cout <<
" -> last cs value in EADL data file (iu) =" << (*pv)[n] <<
G4endl;
250 G4cout <<
"*********************************************************" <<
G4endl;
258 std::vector<G4DynamicParticle*>* fvect,
273 if (verboseLevel > 1) {
274 G4cout <<
"Calling SampleSecondaries() of G4LivermoreNuclearGammaConversionModel"
282 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
285 if (photonEnergy < smallEnergy )
295 if (element ==
nullptr)
297 G4cout <<
"G4LivermoreNuclearGammaConversionModel::SampleSecondaries - element = 0"
302 if (ionisation ==
nullptr)
304 G4cout <<
"G4LivermoreNuclearGammaConversionModel::SampleSecondaries - ionisation = 0"
311 if (photonEnergy > 50. * MeV) fZ += 8. * (element->
GetfCoulomb());
316 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
319 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
320 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
321 G4double epsilonRange = 0.5 - epsilonMin ;
327 G4double f10 = ScreenFunction1(screenMin) - fZ;
328 G4double f20 = ScreenFunction2(screenMin) - fZ;
329 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
330 G4double normF2 = std::max(1.5 * f20,0.);
338 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
344 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
355 electronTotEnergy = (1. -
epsilon) * photonEnergy;
356 positronTotEnergy =
epsilon * photonEnergy;
360 positronTotEnergy = (1. -
epsilon) * photonEnergy;
361 electronTotEnergy =
epsilon * photonEnergy;
381 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
382 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
385 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
386 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
392 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
395 electronDirection.
rotateUz(photonDirection);
402 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
405 positronDirection.
rotateUz(photonDirection);
412 fvect->push_back(particle1);
413 fvect->push_back(particle2);
424G4LivermoreNuclearGammaConversionModel::ScreenFunction1(
G4double screenVariable)
430 if (screenVariable > 1.)
431 value = 42.24 - 8.368 *
G4Log(screenVariable + 0.952);
433 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
441G4LivermoreNuclearGammaConversionModel::ScreenFunction2(
G4double screenVariable)
446 if (screenVariable > 1.)
447 value = 42.24 - 8.368 *
G4Log(screenVariable + 0.952);
449 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
460 G4AutoLock l(&LivermoreNuclearGammaConversionModelMutex);
461 if(!data[
Z]) { ReadData(
Z); }
G4double epsilon(G4double density, G4double temperature)
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
#define G4MUTEX_INITIALIZER
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
G4double GetfCoulomb() const
G4IonisParamElm * GetIonisation() const
G4double GetlogZ3() const
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual ~G4LivermoreNuclearGammaConversionModel()
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4LivermoreNuclearGammaConversionModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LivermoreNuclearConversion")
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
size_t GetNumberOfElements() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
static G4Positron * Positron()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4double HighEnergyLimit() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)