// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4NuclNuclDiffuseElastic.hh 68039 2013-03-13 14:17:40Z gcosmo $ // // // G4 Model: optical elastic scattering with 4-momentum balance // // Class Description // Final state production model for nucleus-nucleus elastic scattering; // Coulomb amplitude is not considered as correction // (as in G4DiffuseElastic) // Class Description - End // // // 17.03.09 V. Grichine implementation for Coulomb elastic scattering #ifndef G4NuclNuclDiffuseElastic_h #define G4NuclNuclDiffuseElastic_h 1 #include #include #include "globals.hh" #include "G4Integrator.hh" #include "G4HadronElastic.hh" #include "G4HadProjectile.hh" #include "G4Nucleus.hh" using namespace std; class G4ParticleDefinition; class G4PhysicsTable; class G4PhysicsLogVector; class G4NuclNuclDiffuseElastic : public G4HadronElastic // G4HadronicInteraction { public: G4NuclNuclDiffuseElastic(); // G4NuclNuclDiffuseElastic(const G4ParticleDefinition* aParticle); virtual ~G4NuclNuclDiffuseElastic(); void Initialise(); void InitialiseOnFly(G4double Z, G4double A); void BuildAngleTable(); virtual G4double SampleInvariantT(const G4ParticleDefinition* p, G4double plab, G4int Z, G4int A); void SetPlabLowLimit(G4double value); void SetHEModelLowLimit(G4double value); void SetQModelLowLimit(G4double value); void SetLowestEnergyLimit(G4double value); void SetRecoilKinEnergyLimit(G4double value); G4double SampleT(const G4ParticleDefinition* aParticle, G4double p, G4double A); G4double SampleTableT(const G4ParticleDefinition* aParticle, G4double p, G4double Z, G4double A); G4double SampleThetaCMS(const G4ParticleDefinition* aParticle, G4double p, G4double A); G4double SampleTableThetaCMS(const G4ParticleDefinition* aParticle, G4double p, G4double Z, G4double A); G4double GetScatteringAngle(G4int iMomentum, G4int iAngle, G4double position); G4double SampleThetaLab(const G4HadProjectile* aParticle, G4double tmass, G4double A); G4double GetDiffuseElasticXsc( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double A ); G4double GetInvElasticXsc( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double A, G4double Z ); G4double GetDiffuseElasticSumXsc( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double A, G4double Z ); G4double GetInvElasticSumXsc( const G4ParticleDefinition* particle, G4double tMand, G4double momentum, G4double A, G4double Z ); G4double IntegralElasticProb( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double A ); G4double GetCoulombElasticXsc( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double Z ); G4double GetRutherfordXsc( G4double theta ); G4double GetInvCoulombElasticXsc( const G4ParticleDefinition* particle, G4double tMand, G4double momentum, G4double A, G4double Z ); G4double GetCoulombTotalXsc( const G4ParticleDefinition* particle, G4double momentum, G4double Z ); G4double GetCoulombIntegralXsc( const G4ParticleDefinition* particle, G4double momentum, G4double Z, G4double theta1, G4double theta2 ); G4double CalculateParticleBeta( const G4ParticleDefinition* particle, G4double momentum ); G4double CalculateZommerfeld( G4double beta, G4double Z1, G4double Z2 ); G4double CalculateAm( G4double momentum, G4double n, G4double Z); G4double CalculateNuclearRad( G4double A); G4double ThetaCMStoThetaLab(const G4DynamicParticle* aParticle, G4double tmass, G4double thetaCMS); G4double ThetaLabToThetaCMS(const G4DynamicParticle* aParticle, G4double tmass, G4double thetaLab); void TestAngleTable(const G4ParticleDefinition* theParticle, G4double partMom, G4double Z, G4double A); G4double BesselJzero(G4double z); G4double BesselJone(G4double z); G4double DampFactor(G4double z); G4double BesselOneByArg(G4double z); G4double GetDiffElasticProb(G4double theta); G4double GetDiffElasticSumProb(G4double theta); G4double GetDiffElasticSumProbA(G4double alpha); G4double GetIntegrandFunction(G4double theta); G4double GetNuclearRadius(){return fNuclearRadius;}; // Technical math functions for strong Coulomb contribution G4complex GammaLogarithm(G4complex xx); G4complex GammaLogB2n(G4complex xx); G4double GetErf(G4double x); G4double GetCosHaPit2(G4double t){return std::cos(CLHEP::halfpi*t*t);}; G4double GetSinHaPit2(G4double t){return std::sin(CLHEP::halfpi*t*t);}; G4double GetCint(G4double x); G4double GetSint(G4double x); G4complex GetErfcComp(G4complex z, G4int nMax); G4complex GetErfcSer(G4complex z, G4int nMax); G4complex GetErfcInt(G4complex z); // , G4int nMax); G4complex GetErfComp(G4complex z, G4int nMax); // AandS algorithm != Ser, Int G4complex GetErfSer(G4complex z, G4int nMax); G4double GetExpCos(G4double x); G4double GetExpSin(G4double x); G4complex GetErfInt(G4complex z); // , G4int nMax); G4double GetLegendrePol(G4int n, G4double x); G4complex TestErfcComp(G4complex z, G4int nMax); G4complex TestErfcSer(G4complex z, G4int nMax); G4complex TestErfcInt(G4complex z); // , G4int nMax); G4complex CoulombAmplitude(G4double theta); G4double CoulombAmplitudeMod2(G4double theta); void CalculateCoulombPhaseZero(); G4double CalculateCoulombPhase(G4int n); void CalculateRutherfordAnglePar(); G4double ProfileNear(G4double theta); G4double ProfileFar(G4double theta); G4double Profile(G4double theta); G4complex PhaseNear(G4double theta); G4complex PhaseFar(G4double theta); G4complex GammaLess(G4double theta); G4complex GammaMore(G4double theta); G4complex AmplitudeNear(G4double theta); G4complex AmplitudeFar(G4double theta); G4complex Amplitude(G4double theta); G4double AmplitudeMod2(G4double theta); G4complex AmplitudeSim(G4double theta); G4double AmplitudeSimMod2(G4double theta); G4double GetRatioSim(G4double theta); G4double GetRatioGen(G4double theta); G4double GetFresnelDiffuseXsc(G4double theta); G4double GetFresnelIntegrandXsc(G4double alpha); G4complex AmplitudeGla(G4double theta); G4double AmplitudeGlaMod2(G4double theta); G4complex AmplitudeGG(G4double theta); G4double AmplitudeGGMod2(G4double theta); void InitParameters(const G4ParticleDefinition* theParticle, G4double partMom, G4double Z, G4double A); void InitDynParameters(const G4ParticleDefinition* theParticle, G4double partMom); void InitParametersGla(const G4DynamicParticle* aParticle, G4double partMom, G4double Z, G4double A); G4double GetHadronNucleonXscNS( G4ParticleDefinition* pParticle, G4double pTkin, G4ParticleDefinition* tParticle); G4double CalcMandelstamS( const G4double mp , const G4double mt , const G4double Plab ); G4double GetProfileLambda(){return fProfileLambda;}; inline void SetProfileLambda(G4double pl) {fProfileLambda = pl;}; inline void SetProfileDelta(G4double pd) {fProfileDelta = pd;}; inline void SetProfileAlpha(G4double pa){fProfileAlpha = pa;}; inline void SetCofLambda(G4double pa){fCofLambda = pa;}; inline void SetCofAlpha(G4double pa){fCofAlpha = pa;}; inline void SetCofAlphaMax(G4double pa){fCofAlphaMax = pa;}; inline void SetCofAlphaCoulomb(G4double pa){fCofAlphaCoulomb = pa;}; inline void SetCofDelta(G4double pa){fCofDelta = pa;}; inline void SetCofPhase(G4double pa){fCofPhase = pa;}; inline void SetCofFar(G4double pa){fCofFar = pa;}; inline void SetEtaRatio(G4double pa){fEtaRatio = pa;}; inline void SetMaxL(G4int l){fMaxL = l;}; inline void SetNuclearRadiusCof(G4double r){fNuclearRadiusCof = r;}; inline G4double GetCofAlphaMax(){return fCofAlphaMax;}; inline G4double GetCofAlphaCoulomb(){return fCofAlphaCoulomb;}; private: G4ParticleDefinition* theProton; G4ParticleDefinition* theNeutron; G4ParticleDefinition* theDeuteron; G4ParticleDefinition* theAlpha; const G4ParticleDefinition* thePionPlus; const G4ParticleDefinition* thePionMinus; G4double lowEnergyRecoilLimit; G4double lowEnergyLimitHE; G4double lowEnergyLimitQ; G4double lowestEnergyLimit; G4double plabLowLimit; G4int fEnergyBin; G4int fAngleBin; G4PhysicsLogVector* fEnergyVector; G4PhysicsTable* fAngleTable; std::vector fAngleBank; std::vector fElementNumberVector; std::vector fElementNameVector; const G4ParticleDefinition* fParticle; G4double fWaveVector; G4double fAtomicWeight; G4double fAtomicNumber; G4double fNuclearRadius1; G4double fNuclearRadius2; G4double fNuclearRadius; G4double fNuclearRadiusSquare; G4double fNuclearRadiusCof; G4double fBeta; G4double fZommerfeld; G4double fRutherfordRatio; G4double fAm; G4bool fAddCoulomb; G4double fCoulombPhase0; G4double fHalfRutThetaTg; G4double fHalfRutThetaTg2; G4double fRutherfordTheta; G4double fProfileLambda; G4double fProfileDelta; G4double fProfileAlpha; G4double fCofLambda; G4double fCofAlpha; G4double fCofDelta; G4double fCofPhase; G4double fCofFar; G4double fCofAlphaMax; G4double fCofAlphaCoulomb; G4int fMaxL; G4double fSumSigma; G4double fEtaRatio; G4double fReZ; }; inline void G4NuclNuclDiffuseElastic::SetRecoilKinEnergyLimit(G4double value) { lowEnergyRecoilLimit = value; } inline void G4NuclNuclDiffuseElastic::SetPlabLowLimit(G4double value) { plabLowLimit = value; } inline void G4NuclNuclDiffuseElastic::SetHEModelLowLimit(G4double value) { lowEnergyLimitHE = value; } inline void G4NuclNuclDiffuseElastic::SetQModelLowLimit(G4double value) { lowEnergyLimitQ = value; } inline void G4NuclNuclDiffuseElastic::SetLowestEnergyLimit(G4double value) { lowestEnergyLimit = value; } //////////////////////////////////////////////////////////////////// // // damp factor in diffraction x/sh(x), x was already *pi inline G4double G4NuclNuclDiffuseElastic::DampFactor(G4double x) { G4double df; G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials // x *= pi; if( std::fabs(x) < 0.01 ) { df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4); } else { df = x/std::sinh(x); } return df; } //////////////////////////////////////////////////////////////////// // // return J1(x)/x with special case for small x inline G4double G4NuclNuclDiffuseElastic::BesselOneByArg(G4double x) { G4double x2, result; if( std::fabs(x) < 0.01 ) { x *= 0.5; x2 = x*x; result = 2. - x2 + x2*x2/6.; } else { result = BesselJone(x)/x; } return result; } //////////////////////////////////////////////////////////////////// // // return particle beta inline G4double G4NuclNuclDiffuseElastic::CalculateParticleBeta(const G4ParticleDefinition* particle, G4double momentum) { G4double mass = particle->GetPDGMass(); G4double a = momentum/mass; fBeta = a/std::sqrt(1+a*a); return fBeta; } //////////////////////////////////////////////////////////////////// // // return Zommerfeld parameter for Coulomb scattering inline G4double G4NuclNuclDiffuseElastic::CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2 ) { fZommerfeld = CLHEP::fine_structure_const*Z1*Z2/beta; return fZommerfeld; } //////////////////////////////////////////////////////////////////// // // return Wentzel correction for Coulomb scattering inline G4double G4NuclNuclDiffuseElastic::CalculateAm( G4double momentum, G4double n, G4double Z) { G4double k = momentum/CLHEP::hbarc; G4double ch = 1.13 + 3.76*n*n; G4double zn = 1.77*k*std::pow(Z,-1./3.)*CLHEP::Bohr_radius; G4double zn2 = zn*zn; fAm = ch/zn2; return fAm; } //////////////////////////////////////////////////////////////////// // // calculate nuclear radius for different atomic weights using different approximations inline G4double G4NuclNuclDiffuseElastic::CalculateNuclearRad( G4double A) { G4double r0 = 1.*CLHEP::fermi, radius; // r0 *= 1.12; // r0 *= 1.44; r0 *= fNuclearRadiusCof; /* if( A < 50. ) { if( A > 10. ) r0 = 1.16*( 1 - std::pow(A, -2./3.) )*CLHEP::fermi; // 1.08*fermi; else r0 = 1.1*CLHEP::fermi; radius = r0*std::pow(A, 1./3.); } else { r0 = 1.7*CLHEP::fermi; // 1.7*fermi; radius = r0*std::pow(A, 0.27); // 0.27); } */ radius = r0*std::pow(A, 1./3.); return radius; } //////////////////////////////////////////////////////////////////// // // return Coulomb scattering differential xsc with Wentzel correction. Test function inline G4double G4NuclNuclDiffuseElastic::GetCoulombElasticXsc( const G4ParticleDefinition* particle, G4double theta, G4double momentum, G4double Z ) { G4double sinHalfTheta = std::sin(0.5*theta); G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta; G4double beta = CalculateParticleBeta( particle, momentum); G4double z = particle->GetPDGCharge(); G4double n = CalculateZommerfeld( beta, z, Z ); G4double am = CalculateAm( momentum, n, Z); G4double k = momentum/CLHEP::hbarc; G4double ch = 0.5*n/k; G4double ch2 = ch*ch; G4double xsc = ch2/((sinHalfTheta2+am)*(sinHalfTheta2+am)); return xsc; } //////////////////////////////////////////////////////////////////// // // return Rutherford scattering differential xsc with Wentzel correction. For Sampling. inline G4double G4NuclNuclDiffuseElastic::GetRutherfordXsc(G4double theta) { G4double sinHalfTheta = std::sin(0.5*theta); G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta; G4double ch2 = fRutherfordRatio*fRutherfordRatio; G4double xsc = ch2/(sinHalfTheta2+fAm)/(sinHalfTheta2+fAm); return xsc; } //////////////////////////////////////////////////////////////////// // // return Coulomb scattering total xsc with Wentzel correction inline G4double G4NuclNuclDiffuseElastic::GetCoulombTotalXsc( const G4ParticleDefinition* particle, G4double momentum, G4double Z ) { G4double beta = CalculateParticleBeta( particle, momentum); G4cout<<"beta = "<GetPDGCharge(); G4double n = CalculateZommerfeld( beta, z, Z ); G4cout<<"fZomerfeld = "<