43 G4int G4BoldyshevTripletModel::maxZ = 99;
61 G4cout <<
"G4BoldyshevTripletModel is constructed " <<
G4endl;
70 for(
G4int i=0; i<maxZ; ++i) {
87 G4cout <<
"Calling Initialise() of G4BoldyshevTripletModel."
104 char* path = getenv(
"G4LEDATA");
111 for(
G4int i=0; i<numOfCouples; ++i)
118 for (
G4int j=0; j<nelm; ++j)
122 else if(Z > maxZ) { Z = maxZ; }
123 if(!
data[Z]) { ReadData(Z, path); }
127 if(isInitialised) {
return; }
129 isInitialised =
true;
147 return lowEnergyLimit;
152 void G4BoldyshevTripletModel::ReadData(
size_t Z,
const char* path)
154 if (verboseLevel > 1)
156 G4cout <<
"Calling ReadData() of G4BoldyshevTripletModel"
160 if(
data[Z]) {
return; }
162 const char* datadir = path;
166 datadir = getenv(
"G4LEDATA");
171 "Environment variable G4LEDATA not defined");
182 std::ostringstream ost;
183 ost << datadir <<
"livermore/tripdata/pp-trip-cs-" << Z <<
".dat";
184 std::ifstream fin(ost.str().c_str());
189 ed <<
"G4BoldyshevTripletModel data file <" << ost.str().c_str()
190 <<
"> is not opened!" <<
G4endl;
193 ed,
"G4LEDATA version should be G4EMLOW6.27 or later.");
200 if(verboseLevel > 3) {
G4cout <<
"File " << ost.str()
201 <<
" is opened by G4BoldyshevTripletModel" <<
G4endl;}
203 data[
Z]->Retrieve(fin,
true);
207 data[
Z] ->SetSpline(
true);
219 if (verboseLevel > 1)
221 G4cout <<
"Calling ComputeCrossSectionPerAtom() of G4BoldyshevTripletModel"
225 if (GammaEnergy < lowEnergyLimit) {
return 0.0; }
231 if(intZ < 1 || intZ > maxZ) {
return xs; }
241 if(!pv) {
return xs; }
244 xs = pv->
Value(GammaEnergy);
249 G4cout <<
"****** DEBUG: tcs value for Z=" << Z <<
" at energy (MeV)="
251 G4cout <<
" cs (Geant4 internal unit)=" << xs <<
G4endl;
252 G4cout <<
" -> first cs value in EADL data file (iu) =" << (*pv)[0] <<
G4endl;
253 G4cout <<
" -> last cs value in EADL data file (iu) =" << (*pv)[
n] <<
G4endl;
254 G4cout <<
"*********************************************************" <<
G4endl;
264 std::vector<G4DynamicParticle*>* fvect,
272 if (verboseLevel > 1) {
273 G4cout <<
"Calling SampleSecondaries() of G4BoldyshevTripletModel"
285 G4double positronTotEnergy, electronTotEnergy, thetaEle, thetaPos;
286 G4double ener_re=0., theta_re, phi_re, phi;
302 G4double SigmaTot = (28./9.) * std::log ( 2.* photonEnergy / electron_mass_c2 ) - 218. / 27. ;
303 G4double X_0 = 2. * ( sqrt(momentumThreshold_N*momentumThreshold_N + 1) -1 );
304 G4double SigmaQ = (82./27. - (14./9.) * log (X_0) + 4./15.*X_0 - 0.0348 * X_0 * X_0);
309 if (recoilProb >= SigmaQ/SigmaTot)
313 ( energyThreshold +
electron_mass_c2 ) / (photonEnergy*momentumThreshold_c) );
320 G4double are, bre, loga, f1_re, greject, cost;
327 cost = pow(cosThetaMax,r1);
328 theta_re = acos(cost);
329 are = 1./(14.*cost*cost);
330 bre = (1.-5.*cost*cost)/(2.*cost);
331 loga = log((1.+ cost)/(1.- cost));
332 f1_re = 1. - bre*loga;
340 }
while(greject < r2);
351 G4double fp = 1. - sint2*loga/(2.*cost) ;
352 rt = (1.-cos(2.*phi_re)*fp/f1_re)/(2.*
pi) ;
360 G4double D2 = 4.*S * electron_mass_c2*electron_mass_c2
362 *(S - electron_mass_c2*electron_mass_c2)*sin(theta_re)*sin(theta_re);
363 ener_re = electron_mass_c2 * (S + electron_mass_c2*
electron_mass_c2)/sqrt(D2);
368 G4double ener_recoil = sqrt( momentum_recoil*momentum_recoil + electron_mass_c2*electron_mass_c2);
369 ener_re = ener_recoil;
374 G4double dxEle_re=sin(theta_re)*std::cos(phi_re),dyEle_re=sin(theta_re)*std::sin(phi_re), dzEle_re=cos(theta_re);
378 G4ThreeVector electronRDirection (dxEle_re, dyEle_re, dzEle_re);
379 electronRDirection.
rotateUz(photonDirection);
383 electronRKineEnergy);
384 fvect->push_back(particle3);
399 G4double t = 0.5*log(momentumThreshold_N + sqrt(momentumThreshold_N*momentumThreshold_N+1));
401 G4double J1 = 0.5*(t*cosh(t)/sinh(t) - log(2.*sinh(t)));
402 G4double J2 = (-2./3.)*log(2.*sinh(t)) + t*cosh(t)/sinh(t) + (sinh(t)-t*pow(cosh(t),3))/(3.*pow(sinh(t),3));
410 G4double b1 = 16. - 3.*b - 36.*b*re*n + 36.*b*pow(re,2.)*pow(n,2.) +
413 G4double c1 = (-6. + 12.*re*n + b + 2*
a)*pow(b,2.);
414 epsilon = (pow(c1,1./3.))/(2.*b) + (b-4.)/(2.*pow(c1,1./3.))+0.5;
416 G4double photonEnergy1 = photonEnergy - ener_re ;
417 positronTotEnergy = epsilon*photonEnergy1;
418 electronTotEnergy = photonEnergy1 - positronTotEnergy;
420 G4double momento_e = sqrt(electronTotEnergy*electronTotEnergy -
421 electron_mass_c2*electron_mass_c2) ;
422 G4double momento_p = sqrt(positronTotEnergy*positronTotEnergy -
423 electron_mass_c2*electron_mass_c2) ;
425 thetaEle = acos((sqrt(p0*p0/(momento_e*momento_e) +1.)- p0/momento_e)) ;
426 thetaPos = acos((sqrt(p0*p0/(momento_p*momento_p) +1.)- p0/momento_p)) ;
429 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
430 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
438 G4double electronKineEnergy =
std::max(0.,electronTotEnergy - electron_mass_c2) ;
442 electronDirection.
rotateUz(photonDirection);
450 G4double positronKineEnergy =
std::max(0.,positronTotEnergy - electron_mass_c2) ;
453 positronDirection.
rotateUz(photonDirection);
458 positronDirection, positronKineEnergy);
461 fvect->push_back(particle1);
462 fvect->push_back(particle2);
485 if(!
data[Z]) { ReadData(Z); }
G4double LowEnergyLimit() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4BoldyshevTripletModel(const G4ParticleDefinition *p=0, const G4String &nam="BoldyshevTripletConversion")
std::vector< G4Element * > G4ElementVector
std::ostringstream G4ExceptionDescription
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4double HighEnergyLimit() const
std::vector< ExP01TrackerHit * > a
virtual ~G4BoldyshevTripletModel()
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
size_t GetVectorLength() const
const G4ElementVector * GetElementVector() const
#define G4MUTEX_INITIALIZER
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
const XML_Char const XML_Char * data
G4GLOB_DLL std::ostream G4cout
size_t GetTableSize() const
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double GeV
static G4Electron * Electron()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static constexpr double MeV
static constexpr double pi
size_t GetNumberOfElements() const
void ProposeTrackStatus(G4TrackStatus status)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
double epsilon(double density, double temperature)
static constexpr double pi
G4ParticleChangeForGamma * GetParticleChangeForGamma()
const G4Material * GetMaterial() const