46 isInitialised(false),meanFreePathTable(0),crossSectionHandler(0)
48 lowEnergyLimit = 2*electron_mass_c2;
49 highEnergyLimit = 100 * GeV;
65 if(verboseLevel > 0) {
66 G4cout <<
"Livermore Polarized GammaConversion is constructed " <<
G4endl
68 << lowEnergyLimit / keV <<
" keV - "
69 << highEnergyLimit / GeV <<
" GeV"
74 crossSectionHandler->
Initialise(0,lowEnergyLimit,highEnergyLimit,400);
81 delete crossSectionHandler;
92 G4cout <<
"Calling G4LivermorePolarizedGammaConversionModel::Initialise()" <<
G4endl;
94 if (crossSectionHandler)
96 crossSectionHandler->
Clear();
97 delete crossSectionHandler;
112 G4cout <<
"G4LivermorePolarizedGammaConversionModel: high energy limit decreased from " <<
121 crossSectionHandler->
Clear();
122 G4String crossSectionFile =
"pair/pp-cs-";
123 crossSectionHandler->
LoadData(crossSectionFile);
126 if (verboseLevel > 2) {
127 G4cout <<
"Loaded cross section files for Livermore Polarized GammaConversion model"
132 if(verboseLevel > 0) {
133 G4cout <<
"Livermore Polarized GammaConversion model is initialized " <<
G4endl
142 isInitialised =
true;
155 if (verboseLevel > 3) {
156 G4cout <<
"G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom()"
159 if(Z < 0.9 || GammaEnergy <= lowEnergyLimit) {
return 0.0; }
179 if (verboseLevel > 3)
180 G4cout <<
"Calling SampleSecondaries() of G4LivermorePolarizedGammaConversionModel" <<
G4endl;
185 if(photonEnergy <= lowEnergyLimit)
199 if(!(gammaPolarization0.
isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.
mag()==0))
201 gammaPolarization0 = GetRandomPolarization(gammaDirection0);
207 gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
215 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
219 if (photonEnergy < smallEnergy )
221 epsilon = epsilon0Local + (0.5 - epsilon0Local) *
G4UniformRand();
254 if (photonEnergy > 50. * MeV) fZ += 8. * (element->
GetfCoulomb());
258 G4double screenMax = exp ((42.24 - fZ)/8.368) - 0.952 ;
259 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
262 G4double epsilon1 = 0.5 - 0.5 * sqrt(1. - screenMin / screenMax) ;
263 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
264 G4double epsilonRange = 0.5 - epsilonMin ;
270 G4double f10 = ScreenFunction1(screenMin) - fZ;
271 G4double f20 = ScreenFunction2(screenMin) - fZ;
272 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
273 G4double normF2 = std::max(1.5 * f20,0.);
278 epsilon = 0.5 - epsilonRange * pow(
G4UniformRand(), 0.3333) ;
279 screen = screenFactor / (epsilon * (1. - epsilon));
280 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
285 screen = screenFactor / (epsilon * (1 - epsilon));
286 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
302 electronTotEnergy = (1. - epsilon) * photonEnergy;
303 positronTotEnergy = epsilon * photonEnergy;
307 positronTotEnergy = (1. - epsilon) * photonEnergy;
308 electronTotEnergy = epsilon * photonEnergy;
330 G4double Ene = electronTotEnergy/electron_mass_c2;
335 SetTheta(&cosTheta,&sinTheta,Ene);
346 phi = SetPhi(photonEnergy);
347 psi = SetPsi(photonEnergy,phi);
372 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
374 SystemOfRefChange(gammaDirection0,electronDirection,
383 Ene = positronTotEnergy/electron_mass_c2;
388 SetTheta(&cosTheta,&sinTheta,Ene);
391 dirX = sinTheta*cos(phip);
392 dirY = sinTheta*sin(phip);
396 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
397 SystemOfRefChange(gammaDirection0,positronDirection,
402 positronDirection, positronKineEnergy);
405 fvect->push_back(particle1);
406 fvect->push_back(particle2);
427G4double G4LivermorePolarizedGammaConversionModel::ScreenFunction1(
G4double screenVariable)
433 if (screenVariable > 1.)
434 value = 42.24 - 8.368 * log(screenVariable + 0.952);
436 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
443G4double G4LivermorePolarizedGammaConversionModel::ScreenFunction2(
G4double screenVariable)
449 if (screenVariable > 1.)
450 value = 42.24 - 8.368 * log(screenVariable + 0.952);
452 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
464 G4double Momentum = sqrt(Energy*Energy -1);
467 *p_cosTheta = (Energy*((2*Rand)- 1) + Momentum)/((Momentum*(2*Rand-1))+Energy);
468 *p_sinTheta = (2*sqrt(Rand*(1-Rand)))/(Momentum*(2*Rand-1)+Energy);
492 const G4double ay0=5.6, by0=18.6, aa0=2.9, ba0 = 8.16E-3;
493 const G4double aw = 0.0151, bw = 10.7, cw = -410.;
495 const G4double axc = 3.1455, bxc = -1.11, cxc = 310.;
497 pl[0] = Fln(ay0,by0,Ene);
498 pl[1] = aa0 + ba0*(Ene);
499 pl[2] = Poli(aw,bw,cw,Ene);
500 pl[3] = Poli(axc,bxc,cxc,Ene);
502 const G4double abf = 3.1216, bbf = 2.68;
504 pt[1] = abf + bbf/Ene;
513 n1 = Fintlor(pl,pi) - Fintlor(pl,xe);
514 n2 = Finttan(pt,xe) - Finttan(pt,0.);
520 const G4double aw = 0.21, bw = 10.8, cw = -58.;
521 const G4double axc = 3.17, bxc = -0.87, cxc = -6.;
523 pl[0] = Fln(ay0, by0, Ene);
524 pl[1] = Fln(aa0, ba0, Ene);
525 pl[2] = Poli(aw,bw,cw,Ene);
526 pl[3] = Poli(axc,bxc,cxc,Ene);
530 n1 = Fintlor(pl,pi) - Fintlor(pl,xe);
557 value = Finvlor(pl,xe,r2);
558 xco = Glor(pl,value)/c1;
564 value = Finvtan(pt,n,r1);
572 value = Finvlor(pl,xe,r2);
573 xco = Glor(pl,value)/c1;
598 const G4double ay00 = 3.4, by00 = 9.8, aa00 = 1.34, ba00 = 5.3;
599 const G4double aw0 = 0.014, bw0 = 9.7, cw0 = -2.E4;
600 const G4double axc0 = 3.1423, bxc0 = -2.35, cxc0 = 0.;
601 const G4double ay0p = 1.53, by0p = 3.2, aap = 0.67, bap = 8.5E-3;
602 const G4double awp = 6.9E-3, bwp = 12.6, cwp = -3.8E4;
603 const G4double axcp = 2.8E-3,bxcp = -3.133;
604 const G4double abf0 = 3.1213, bbf0 = 2.61;
605 const G4double abfpm = 3.1231, bbfpm = 2.84;
607 p0l[0] = Fln(ay00, by00, Ene);
608 p0l[1] = Fln(aa00, ba00, Ene);
609 p0l[2] = Poli(aw0, bw0, cw0, Ene);
610 p0l[3] = Poli(axc0, bxc0, cxc0, Ene);
612 ppml[0] = Fln(ay0p, by0p, Ene);
613 ppml[1] = aap + bap*(Ene);
614 ppml[2] = Poli(awp, bwp, cwp, Ene);
615 ppml[3] = Fln(axcp,bxcp,Ene);
618 p0t[1] = abf0 + bbf0/Ene;
620 ppmt[1] = abfpm + bbfpm/Ene;
626 xe0 = Encu(p0l, p0t, xi);
628 xepm = Encu(ppml, ppmt, xi);
633 const G4double ay00 = 2.82, by00 = 6.35;
634 const G4double aa00 = -1.75, ba00 = 0.25;
636 const G4double aw0 = 0.028, bw0 = 5., cw0 = -50.;
637 const G4double axc0 = 3.14213, bxc0 = -2.3, cxc0 = 5.7;
638 const G4double ay0p = 1.56, by0p = 3.6;
639 const G4double aap = 0.86, bap = 8.3E-3;
640 const G4double awp = 0.022, bwp = 7.4, cwp = -51.;
643 p0l[0] = Fln(ay00, by00, Ene);
644 p0l[1] = aa00+pow(Ene, ba00);
645 p0l[2] = Poli(aw0, bw0, cw0, Ene);
646 p0l[3] = Poli(axc0, bxc0, cxc0, Ene);
647 ppml[0] = Fln(ay0p, by0p, Ene);
648 ppml[1] = aap + bap*(Ene);
649 ppml[2] = Poli(awp, bwp, cwp, Ene);
660 b = (ppml[0]+2*ppml[1]*ppml[2]*Flor(ppml,PhiLocal));
664 b = Ftan(ppmt,PhiLocal);
668 a = (p0l[0]+2*p0l[1]*p0l[2]*Flor(p0l,PhiLocal));
672 a = Ftan(p0t,PhiLocal);
677 b = (ppml[0]+2*ppml[1]*ppml[2]*Flor(ppml,PhiLocal));
678 a = (p0l[0]+2*p0l[1]*p0l[2]*Flor(p0l,PhiLocal));
698 r3 = nr*(a*cos(value)*cos(value) + b*sin(value)*sin(value));
705G4double G4LivermorePolarizedGammaConversionModel::Poli
711 value =(a + b/x + c/(x*x*x));
719G4double G4LivermorePolarizedGammaConversionModel::Fln
735G4double G4LivermorePolarizedGammaConversionModel::Encu
745 fx = (Flor(p_p1,x)*Glor(p_p1,x) - Ftan(p_p2, x))/
746 (Fdlor(p_p1,x) - Fdtan(p_p2,x));
748 if(x > xmax) {
return xmax; }
753 if(std::fabs(fx) <= x*1.0e-6) {
break; }
756 if(x < 0.0) { x = 0.0; }
769 value = 1./(
pi*(w*w + 4.*(x-xc)*(x-xc)));
782 value = (y0 *
pi*(w*w + 4.*(x-xc)*(x-xc)) + 2.*A*w);
795 value = (-16.*A*w*(x-xc))/
796 (
pi*(w*w+4.*(x-xc)*(x-xc))*(w*w+4.*(x-xc)*(x-xc)));
809 value = y0*x + A*atan( 2*(x-xc)/w) /
pi;
823 nor = atan(2.*(pi-xc)/w)/(2.*
pi*w) - atan(2.*(x-xc)/w)/(2.*pi*w);
824 value = xc - (w/2.)*tan(-2.*r*nor*pi*w+atan(2*(xc-x)/w));
847 value = -1.*a / ((x-b)*(x-b));
870 value = b*(1-exp(r*cnor/a));
910 c.
setX(std::cos(angle)*(a0.
x())+std::sin(angle)*b0.
x());
911 c.
setY(std::cos(angle)*(a0.
y())+std::sin(angle)*b0.
y());
912 c.
setZ(std::cos(angle)*(a0.
z())+std::sin(angle)*b0.
z());
922G4ThreeVector G4LivermorePolarizedGammaConversionModel::GetPerpendicularPolarization
936 return gammaPolarization - gammaPolarization.
dot(gammaDirection)/gammaDirection.
dot(gammaDirection) * gammaDirection;
942void G4LivermorePolarizedGammaConversionModel::SystemOfRefChange
957 direction1 = (direction_x*Axis_X0 + direction_y*Axis_Y0 + direction_z*Axis_Z0).unit();
CLHEP::Hep3Vector G4ThreeVector
G4DLLIMPORT std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
double howOrthogonal(const Hep3Vector &v) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static G4Electron * Electron()
G4double GetfCoulomb() const
G4IonisParamElm * GetIonisation() const
G4double GetlogZ3() const
G4LivermorePolarizedGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedGammaConversion")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChange
virtual ~G4LivermorePolarizedGammaConversionModel()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static G4Positron * Positron()
G4double FindValue(G4int Z, G4double e) const
void LoadData(const G4String &dataFile)
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
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 SetLowEnergyLimit(G4double)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)