Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DiffuseElastic Class Reference

#include <G4DiffuseElastic.hh>

Inheritance diagram for G4DiffuseElastic:

Public Member Functions

 G4DiffuseElastic ()
virtual ~G4DiffuseElastic ()
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
void Initialise ()
void InitialiseOnFly (G4double Z, G4double A)
void BuildAngleTable ()
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4double NeutronTuniform (G4int Z)
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 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 ()
Public Member Functions inherited from G4HadronElastic
 G4HadronElastic (const G4String &name="hElasticLHEP")
 ~G4HadronElastic () override
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
G4double GetSlopeCof (const G4int pdg)
void SetLowestEnergyLimit (G4double value)
G4double LowestEnergyLimit () const
G4double ComputeMomentumCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
void ModelDescription (std::ostream &) const override
Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
virtual ~G4HadronicInteraction ()
G4double GetMinEnergy () const
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
void SetMinEnergy (G4double anEnergy)
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
G4double GetMaxEnergy () const
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
void SetMaxEnergy (const G4double anEnergy)
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
G4int GetVerboseLevel () const
void SetVerboseLevel (G4int value)
const G4StringGetModelName () const
void DeActivateFor (const G4Material *aMaterial)
void ActivateFor (const G4Material *aMaterial)
void DeActivateFor (const G4Element *anElement)
void ActivateFor (const G4Element *anElement)
G4bool IsBlocked (const G4Material *aMaterial) const
G4bool IsBlocked (const G4Element *anElement) const
void SetRecoilEnergyThreshold (G4double val)
G4double GetRecoilEnergyThreshold () const
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual void InitialiseModel ()
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
G4bool operator== (const G4HadronicInteraction &right) const =delete
G4bool operator!= (const G4HadronicInteraction &right) const =delete

Additional Inherited Members

Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
G4bool IsBlocked () const
void Block ()
Protected Attributes inherited from G4HadronElastic
G4double pLocalTmax
G4int secID
Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
G4int verboseLevel
G4double theMinEnergy
G4double theMaxEnergy
G4bool isBlocked

Detailed Description

Definition at line 58 of file G4DiffuseElastic.hh.

Constructor & Destructor Documentation

◆ G4DiffuseElastic()

G4DiffuseElastic::G4DiffuseElastic ( )

Definition at line 75 of file G4DiffuseElastic.cc.

76 : G4HadronElastic("DiffuseElastic"), fParticle(0)
77{
78 SetMinEnergy( 0.01*MeV );
80
81 verboseLevel = 0;
82 lowEnergyRecoilLimit = 100.*keV;
83 lowEnergyLimitQ = 0.0*GeV;
84 lowEnergyLimitHE = 0.0*GeV;
85 lowestEnergyLimit = 0.0*keV;
86 plabLowLimit = 20.0*MeV;
87
88 theProton = G4Proton::Proton();
89 theNeutron = G4Neutron::Neutron();
90 theDeuteron = G4Deuteron::Deuteron();
91 theAlpha = G4Alpha::Alpha();
92 thePionPlus = G4PionPlus::PionPlus();
93 thePionMinus = G4PionMinus::PionMinus();
94
95 fEnergyBin = 300; // Increased from the original 200 to have no wider log-energy-bins up to 10 PeV
96 fAngleBin = 200;
97
98 fEnergyVector = new G4PhysicsLogVector( theMinEnergy, theMaxEnergy, fEnergyBin );
99
100 fAngleTable = 0;
101
102 fParticle = 0;
103 fWaveVector = 0.;
104 fAtomicWeight = 0.;
105 fAtomicNumber = 0.;
106 fNuclearRadius = 0.;
107 fBeta = 0.;
108 fZommerfeld = 0.;
109 fAm = 0.;
110 fAddCoulomb = false;
111}
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
G4HadronElastic(const G4String &name="hElasticLHEP")
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4Proton * Proton()
Definition G4Proton.cc:90

Referenced by BuildAngleTable(), IntegralElasticProb(), SampleThetaCMS(), and TestAngleTable().

◆ ~G4DiffuseElastic()

G4DiffuseElastic::~G4DiffuseElastic ( )
virtual

Definition at line 117 of file G4DiffuseElastic.cc.

118{
119 if ( fEnergyVector )
120 {
121 delete fEnergyVector;
122 fEnergyVector = 0;
123 }
124 for ( std::vector<G4PhysicsTable*>::iterator it = fAngleBank.begin();
125 it != fAngleBank.end(); ++it )
126 {
127 if ( (*it) ) (*it)->clearAndDestroy();
128
129 delete *it;
130 *it = 0;
131 }
132 fAngleTable = 0;
133}

Member Function Documentation

◆ BesselJone()

G4double G4DiffuseElastic::BesselJone ( G4double z)
inline

Definition at line 329 of file G4DiffuseElastic.hh.

330{
331 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
332
333 modvalue = std::fabs(value);
334
335 if ( modvalue < 8.0 )
336 {
337 value2 = value*value;
338
339 fact1 = value*(72362614232.0 + value2*(-7895059235.0
340 + value2*( 242396853.1
341 + value2*(-2972611.439
342 + value2*( 15704.48260
343 + value2*(-30.16036606 ) ) ) ) ) );
344
345 fact2 = 144725228442.0 + value2*(2300535178.0
346 + value2*(18583304.74
347 + value2*(99447.43394
348 + value2*(376.9991397
349 + value2*1.0 ) ) ) );
350 bessel = fact1/fact2;
351 }
352 else
353 {
354 arg = 8.0/modvalue;
355
356 value2 = arg*arg;
357
358 shift = modvalue - 2.356194491;
359
360 fact1 = 1.0 + value2*( 0.183105e-2
361 + value2*(-0.3516396496e-4
362 + value2*(0.2457520174e-5
363 + value2*(-0.240337019e-6 ) ) ) );
364
365 fact2 = 0.04687499995 + value2*(-0.2002690873e-3
366 + value2*( 0.8449199096e-5
367 + value2*(-0.88228987e-6
368 + value2*0.105787412e-6 ) ) );
369
370 bessel = std::sqrt( 0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
371
372 if (value < 0.0) bessel = -bessel;
373 }
374 return bessel;
375}
double G4double
Definition G4Types.hh:83

Referenced by BesselOneByArg(), GetDiffElasticProb(), GetDiffElasticSumProb(), and GetDiffElasticSumProbA().

◆ BesselJzero()

G4double G4DiffuseElastic::BesselJzero ( G4double z)
inline

Definition at line 277 of file G4DiffuseElastic.hh.

278{
279 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
280
281 modvalue = std::fabs(value);
282
283 if ( value < 8.0 && value > -8.0 )
284 {
285 value2 = value*value;
286
287 fact1 = 57568490574.0 + value2*(-13362590354.0
288 + value2*( 651619640.7
289 + value2*(-11214424.18
290 + value2*( 77392.33017
291 + value2*(-184.9052456 ) ) ) ) );
292
293 fact2 = 57568490411.0 + value2*( 1029532985.0
294 + value2*( 9494680.718
295 + value2*(59272.64853
296 + value2*(267.8532712
297 + value2*1.0 ) ) ) );
298
299 bessel = fact1/fact2;
300 }
301 else
302 {
303 arg = 8.0/modvalue;
304
305 value2 = arg*arg;
306
307 shift = modvalue-0.785398164;
308
309 fact1 = 1.0 + value2*(-0.1098628627e-2
310 + value2*(0.2734510407e-4
311 + value2*(-0.2073370639e-5
312 + value2*0.2093887211e-6 ) ) );
313
314 fact2 = -0.1562499995e-1 + value2*(0.1430488765e-3
315 + value2*(-0.6911147651e-5
316 + value2*(0.7621095161e-6
317 - value2*0.934945152e-7 ) ) );
318
319 bessel = std::sqrt(0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2 );
320 }
321 return bessel;
322}

Referenced by GetDiffElasticProb(), GetDiffElasticSumProb(), and GetDiffElasticSumProbA().

◆ BesselOneByArg()

G4double G4DiffuseElastic::BesselOneByArg ( G4double z)
inline

Definition at line 404 of file G4DiffuseElastic.hh.

405{
406 G4double x2, result;
407
408 if( std::fabs(x) < 0.01 )
409 {
410 x *= 0.5;
411 x2 = x*x;
412 result = 2. - x2 + x2*x2/6.;
413 }
414 else
415 {
416 result = BesselJone(x)/x;
417 }
418 return result;
419}
G4double BesselJone(G4double z)

Referenced by GetDiffElasticProb(), GetDiffElasticSumProb(), and GetDiffElasticSumProbA().

◆ BuildAngleTable()

void G4DiffuseElastic::BuildAngleTable ( )

Definition at line 1001 of file G4DiffuseElastic.cc.

1002{
1003 G4int i, j;
1004 G4double partMom, kinE, a = 0., z = fParticle->GetPDGCharge(), m1 = fParticle->GetPDGMass();
1005 G4double alpha1, alpha2, alphaMax, alphaCoulomb, delta = 0., sum = 0.;
1006
1007 G4Integrator<G4DiffuseElastic,G4double(G4DiffuseElastic::*)(G4double)> integral;
1008
1009 fAngleTable = new G4PhysicsTable( fEnergyBin );
1010
1011 for( i = 0; i < fEnergyBin; i++)
1012 {
1013 kinE = fEnergyVector->GetLowEdgeEnergy(i);
1014 partMom = std::sqrt( kinE*(kinE + 2*m1) );
1015
1016 fWaveVector = partMom/hbarc;
1017
1018 G4double kR = fWaveVector*fNuclearRadius;
1019 G4double kR2 = kR*kR;
1020 G4double kRmax = 18.6; // 10.6; 10.6, 18, 10.174; ~ 3 maxima of J1 or 15., 25.
1021 G4double kRcoul = 1.9; // 1.2; 1.4, 2.5; // on the first slope of J1
1022 // G4double kRlim = 1.2;
1023 // G4double kRlim2 = kRlim*kRlim/kR2;
1024
1025 alphaMax = kRmax*kRmax/kR2;
1026
1027
1028 // if (alphaMax > 4.) alphaMax = 4.; // vmg05-02-09: was pi2
1029 // if ( alphaMax > 4. || alphaMax < 1. ) alphaMax = 15.; // vmg27.11.14
1030
1031 // if ( alphaMax > 4. || alphaMax < 1. ) alphaMax = CLHEP::pi*CLHEP::pi; // vmg06.01.15
1032
1033 // G4cout<<"alphaMax = "<<alphaMax<<", ";
1034
1035 if ( alphaMax >= CLHEP::pi*CLHEP::pi ) alphaMax = CLHEP::pi*CLHEP::pi; // vmg21.10.15
1036
1037 alphaCoulomb = kRcoul*kRcoul/kR2;
1038
1039 if( z )
1040 {
1041 a = partMom/m1; // beta*gamma for m1
1042 fBeta = a/std::sqrt(1+a*a);
1043 fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
1044 fAm = CalculateAm( partMom, fZommerfeld, fAtomicNumber);
1045 }
1046 G4PhysicsFreeVector* angleVector = new G4PhysicsFreeVector(fAngleBin-1);
1047
1048 // G4PhysicsLogVector* angleBins = new G4PhysicsLogVector( 0.001*alphaMax, alphaMax, fAngleBin );
1049
1050 G4double delth = alphaMax/fAngleBin;
1051
1052 sum = 0.;
1053
1054 // fAddCoulomb = false;
1055 fAddCoulomb = true;
1056
1057 // for(j = 1; j < fAngleBin; j++)
1058 for(j = fAngleBin-1; j >= 1; j--)
1059 {
1060 // alpha1 = angleBins->GetLowEdgeEnergy(j-1);
1061 // alpha2 = angleBins->GetLowEdgeEnergy(j);
1062
1063 // alpha1 = alphaMax*(j-1)/fAngleBin;
1064 // alpha2 = alphaMax*( j )/fAngleBin;
1065
1066 alpha1 = delth*(j-1);
1067 // if(alpha1 < kRlim2) alpha1 = kRlim2;
1068 alpha2 = alpha1 + delth;
1069
1070 // if( ( alpha2 > alphaCoulomb ) && z ) fAddCoulomb = true;
1071 if( ( alpha1 < alphaCoulomb ) && z ) fAddCoulomb = false;
1072
1073 delta = integral.Legendre10(this, &G4DiffuseElastic::GetIntegrandFunction, alpha1, alpha2);
1074 // delta = integral.Legendre96(this, &G4DiffuseElastic::GetIntegrandFunction, alpha1, alpha2);
1075
1076 sum += delta;
1077
1078 angleVector->PutValue( j-1 , alpha1, sum ); // alpha2
1079 // G4cout<<"j-1 = "<<j-1<<"; alpha2 = "<<alpha2 << " delta "<< delta <<"; sum = "<<sum<<G4endl;
1080 }
1081 fAngleTable->insertAt(i, angleVector);
1082
1083 // delete[] angleVector;
1084 // delete[] angleBins;
1085 }
1086 return;
1087}
int G4int
Definition G4Types.hh:85
const G4double alpha2
G4double CalculateAm(G4double momentum, G4double n, G4double Z)
G4double GetIntegrandFunction(G4double theta)
G4double CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2)
void PutValue(const std::size_t index, const G4double e, const G4double value)

Referenced by Initialise(), and InitialiseOnFly().

◆ CalculateAm()

G4double G4DiffuseElastic::CalculateAm ( G4double momentum,
G4double n,
G4double Z )
inline

Definition at line 450 of file G4DiffuseElastic.hh.

451{
452 G4double k = momentum/CLHEP::hbarc;
453 G4double ch = 1.13 + 3.76*n*n;
454 G4double zn = 1.77*k*(1.0/G4Pow::GetInstance()->A13(Z))*CLHEP::Bohr_radius;
455 G4double zn2 = zn*zn;
456 fAm = ch/zn2;
457
458 return fAm;
459}
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double A13(G4double A) const
Definition G4Pow.cc:116

Referenced by BuildAngleTable(), GetCoulombElasticXsc(), GetCoulombIntegralXsc(), GetCoulombTotalXsc(), GetDiffuseElasticSumXsc(), and TestAngleTable().

◆ CalculateNuclearRad()

G4double G4DiffuseElastic::CalculateNuclearRad ( G4double A)
inline

Definition at line 465 of file G4DiffuseElastic.hh.

466{
467 G4double R, r0, a11, a12, a13, a2, a3;
468
469 a11 = 1.26; // 1.08, 1.16
470 a12 = 1.; // 1.08, 1.16
471 a13 = 1.12; // 1.08, 1.16
472 a2 = 1.1;
473 a3 = 1.;
474
475 // Special rms radii for light nucleii
476
477 if (A < 50.)
478 {
479 if (std::abs(A-1.) < 0.5) return 0.89*CLHEP::fermi; // p
480 else if(std::abs(A-2.) < 0.5) return 2.13*CLHEP::fermi; // d
481 else if( // std::abs(Z-1.) < 0.5 &&
482std::abs(A-3.) < 0.5) return 1.80*CLHEP::fermi; // t
483
484 // else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96CLHEP::fermi; // He3
485 else if( // std::abs(Z-2.) < 0.5 &&
486std::abs(A-4.) < 0.5) return 1.68*CLHEP::fermi; // He4
487
488 else if( // std::abs(Z-3.) < 0.5
489 std::abs(A-7.) < 0.5 ) return 2.40*CLHEP::fermi; // Li7
490 else if( // std::abs(Z-4.) < 0.5
491std::abs(A-9.) < 0.5) return 2.51*CLHEP::fermi; // Be9
492
493 else if( 10. < A && A <= 16. ) r0 = a11*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi; // 1.08CLHEP::fermi;
494 else if( 15. < A && A <= 20. ) r0 = a12*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
495 else if( 20. < A && A <= 30. ) r0 = a13*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
496 else r0 = a2*CLHEP::fermi;
497
498 R = r0*G4Pow::GetInstance()->A13(A);
499 }
500 else
501 {
502 r0 = a3*CLHEP::fermi;
503
504 R = r0*G4Pow::GetInstance()->powA(A, 0.27);
505 }
506 fNuclearRadius = R;
507 return R;
508 /*
509 G4double r0;
510 if( A < 50. )
511 {
512 if( A > 10. ) r0 = 1.16*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi; // 1.08CLHEP::fermi;
513 else r0 = 1.1*CLHEP::fermi;
514 fNuclearRadius = r0*G4Pow::GetInstance()->A13(A);
515 }
516 else
517 {
518 r0 = 1.7*CLHEP::fermi; // 1.7*CLHEP::fermi;
519 fNuclearRadius = r0*G4Pow::GetInstance()->powA(A, 0.27); // 0.27);
520 }
521 return fNuclearRadius;
522 */
523}
const G4double A[17]
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
G4double A23(G4double A) const
Definition G4Pow.hh:131

Referenced by GetDiffuseElasticSumXsc(), GetDiffuseElasticXsc(), Initialise(), InitialiseOnFly(), IntegralElasticProb(), SampleThetaCMS(), and TestAngleTable().

◆ CalculateParticleBeta()

G4double G4DiffuseElastic::CalculateParticleBeta ( const G4ParticleDefinition * particle,
G4double momentum )
inline

Definition at line 425 of file G4DiffuseElastic.hh.

427{
428 G4double mass = particle->GetPDGMass();
429 G4double a = momentum/mass;
430 fBeta = a/std::sqrt(1+a*a);
431
432 return fBeta;
433}

Referenced by GetCoulombElasticXsc(), GetCoulombIntegralXsc(), GetCoulombTotalXsc(), and GetDiffuseElasticSumXsc().

◆ CalculateZommerfeld()

G4double G4DiffuseElastic::CalculateZommerfeld ( G4double beta,
G4double Z1,
G4double Z2 )
inline

Definition at line 439 of file G4DiffuseElastic.hh.

440{
441 fZommerfeld = CLHEP::fine_structure_const*Z1*Z2/beta;
442
443 return fZommerfeld;
444}

Referenced by BuildAngleTable(), GetCoulombElasticXsc(), GetCoulombIntegralXsc(), GetCoulombTotalXsc(), GetDiffuseElasticSumXsc(), and TestAngleTable().

◆ DampFactor()

G4double G4DiffuseElastic::DampFactor ( G4double z)
inline

Definition at line 381 of file G4DiffuseElastic.hh.

382{
383 G4double df;
384 G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials
385
386 // x *= pi;
387
388 if( std::fabs(x) < 0.01 )
389 {
390 df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4);
391 }
392 else
393 {
394 df = x/std::sinh(x);
395 }
396 return df;
397}

Referenced by GetDiffElasticProb(), GetDiffElasticSumProb(), and GetDiffElasticSumProbA().

◆ GetCoulombElasticXsc()

G4double G4DiffuseElastic::GetCoulombElasticXsc ( const G4ParticleDefinition * particle,
G4double theta,
G4double momentum,
G4double Z )
inline

Definition at line 529 of file G4DiffuseElastic.hh.

533{
534 G4double sinHalfTheta = std::sin(0.5*theta);
535 G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
536 G4double beta = CalculateParticleBeta( particle, momentum);
537 G4double z = particle->GetPDGCharge();
538 G4double n = CalculateZommerfeld( beta, z, Z );
539 G4double am = CalculateAm( momentum, n, Z);
540 G4double k = momentum/CLHEP::hbarc;
541 G4double ch = 0.5*n/k;
542 G4double ch2 = ch*ch;
543 G4double xsc = ch2/(sinHalfTheta2+am)/(sinHalfTheta2+am);
544
545 return xsc;
546}
G4double CalculateParticleBeta(const G4ParticleDefinition *particle, G4double momentum)

Referenced by GetInvCoulombElasticXsc().

◆ GetCoulombIntegralXsc()

G4double G4DiffuseElastic::GetCoulombIntegralXsc ( const G4ParticleDefinition * particle,
G4double momentum,
G4double Z,
G4double theta1,
G4double theta2 )
inline

Definition at line 578 of file G4DiffuseElastic.hh.

581{
582 G4double c1 = std::cos(theta1);
583 G4cout<<"c1 = "<<c1<<G4endl;
584 G4double c2 = std::cos(theta2);
585 G4cout<<"c2 = "<<c2<<G4endl;
586 G4double beta = CalculateParticleBeta( particle, momentum);
587 // G4cout<<"beta = "<<beta<<G4endl;
588 G4double z = particle->GetPDGCharge();
589 G4double n = CalculateZommerfeld( beta, z, Z );
590 // G4cout<<"fZomerfeld = "<<n<<G4endl;
591 G4double am = CalculateAm( momentum, n, Z);
592 // G4cout<<"cof Am = "<<am<<G4endl;
593 G4double k = momentum/CLHEP::hbarc;
594 // G4cout<<"k = "<<k*CLHEP::fermi<<" 1/fermi"<<G4endl;
595 // G4cout<<"k*Bohr_radius = "<<k*CLHEP::Bohr_radius<<G4endl;
596 G4double ch = n/k;
597 G4double ch2 = ch*ch;
598 am *= 2.;
599 G4double xsc = ch2*CLHEP::twopi*(c1-c2);
600 xsc /= (1 - c1 + am)*(1 - c2 + am);
601
602 return xsc;
603}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ GetCoulombTotalXsc()

G4double G4DiffuseElastic::GetCoulombTotalXsc ( const G4ParticleDefinition * particle,
G4double momentum,
G4double Z )
inline

Definition at line 553 of file G4DiffuseElastic.hh.

555{
556 G4double beta = CalculateParticleBeta( particle, momentum);
557 G4cout<<"beta = "<<beta<<G4endl;
558 G4double z = particle->GetPDGCharge();
559 G4double n = CalculateZommerfeld( beta, z, Z );
560 G4cout<<"fZomerfeld = "<<n<<G4endl;
561 G4double am = CalculateAm( momentum, n, Z);
562 G4cout<<"cof Am = "<<am<<G4endl;
563 G4double k = momentum/CLHEP::hbarc;
564 G4cout<<"k = "<<k*CLHEP::fermi<<" 1/fermi"<<G4endl;
565 G4cout<<"k*Bohr_radius = "<<k*CLHEP::Bohr_radius<<G4endl;
566 G4double ch = n/k;
567 G4double ch2 = ch*ch;
568 G4double xsc = ch2*CLHEP::pi/(am +am*am);
569
570 return xsc;
571}

◆ GetDiffElasticProb()

G4double G4DiffuseElastic::GetDiffElasticProb ( G4double theta)

Definition at line 380 of file G4DiffuseElastic.cc.

385{
386 G4double sigma, bzero, bzero2, bonebyarg, bonebyarg2, damp, damp2;
387 G4double delta, diffuse, gamma;
388 G4double e1, e2, bone, bone2;
389
390 // G4double wavek = momentum/hbarc; // wave vector
391 // G4double r0 = 1.08*fermi;
392 // G4double rad = r0*G4Pow::GetInstance()->A13(A);
393
394 if (fParticle == theProton)
395 {
396 diffuse = 0.63*fermi;
397 gamma = 0.3*fermi;
398 delta = 0.1*fermi*fermi;
399 e1 = 0.3*fermi;
400 e2 = 0.35*fermi;
401 }
402 else if (fParticle == theNeutron)
403 {
404 diffuse = 0.63*fermi; // 1.63*fermi; //
405 G4double k0 = 1*GeV/hbarc;
406 diffuse *= k0/fWaveVector;
407
408 gamma = 0.3*fermi;
409 delta = 0.1*fermi*fermi;
410 e1 = 0.3*fermi;
411 e2 = 0.35*fermi;
412 }
413 else // as proton, if were not defined
414 {
415 diffuse = 0.63*fermi;
416 gamma = 0.3*fermi;
417 delta = 0.1*fermi*fermi;
418 e1 = 0.3*fermi;
419 e2 = 0.35*fermi;
420 }
421 G4double kr = fWaveVector*fNuclearRadius; // wavek*rad;
422 G4double kr2 = kr*kr;
423 G4double krt = kr*theta;
424
425 bzero = BesselJzero(krt);
426 bzero2 = bzero*bzero;
427 bone = BesselJone(krt);
428 bone2 = bone*bone;
429 bonebyarg = BesselOneByArg(krt);
430 bonebyarg2 = bonebyarg*bonebyarg;
431
432 G4double lambda = 15.; // 15 ok
433
434 // G4double kgamma = fWaveVector*gamma; // wavek*delta;
435
436 G4double kgamma = lambda*(1.-G4Exp(-fWaveVector*gamma/lambda)); // wavek*delta;
437 G4double kgamma2 = kgamma*kgamma;
438
439 // G4double dk2t = delta*fWaveVector*fWaveVector*theta; // delta*wavek*wavek*theta;
440 // G4double dk2t2 = dk2t*dk2t;
441 // G4double pikdt = pi*fWaveVector*diffuse*theta;// pi*wavek*diffuse*theta;
442
443 G4double pikdt = lambda*(1.-G4Exp(-pi*fWaveVector*diffuse*theta/lambda)); // wavek*delta;
444
445 damp = DampFactor(pikdt);
446 damp2 = damp*damp;
447
448 G4double mode2k2 = (e1*e1+e2*e2)*fWaveVector*fWaveVector;
449 G4double e2dk3t = -2.*e2*delta*fWaveVector*fWaveVector*fWaveVector*theta;
450
451
452 sigma = kgamma2;
453 // sigma += dk2t2;
454 sigma *= bzero2;
455 sigma += mode2k2*bone2 + e2dk3t*bzero*bone;
456 sigma += kr2*bonebyarg2;
457 sigma *= damp2; // *rad*rad;
458
459 return sigma;
460}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double BesselOneByArg(G4double z)
G4double BesselJzero(G4double z)
G4double DampFactor(G4double z)

Referenced by GetDiffuseElasticXsc().

◆ GetDiffElasticSumProb()

G4double G4DiffuseElastic::GetDiffElasticSumProb ( G4double theta)

Definition at line 468 of file G4DiffuseElastic.cc.

473{
474 G4double sigma, bzero, bzero2, bonebyarg, bonebyarg2, damp, damp2;
475 G4double delta, diffuse, gamma;
476 G4double e1, e2, bone, bone2;
477
478 // G4double wavek = momentum/hbarc; // wave vector
479 // G4double r0 = 1.08*fermi;
480 // G4double rad = r0*G4Pow::GetInstance()->A13(A);
481
482 G4double kr = fWaveVector*fNuclearRadius; // wavek*rad;
483 G4double kr2 = kr*kr;
484 G4double krt = kr*theta;
485
486 bzero = BesselJzero(krt);
487 bzero2 = bzero*bzero;
488 bone = BesselJone(krt);
489 bone2 = bone*bone;
490 bonebyarg = BesselOneByArg(krt);
491 bonebyarg2 = bonebyarg*bonebyarg;
492
493 if (fParticle == theProton)
494 {
495 diffuse = 0.63*fermi;
496 // diffuse = 0.6*fermi;
497 gamma = 0.3*fermi;
498 delta = 0.1*fermi*fermi;
499 e1 = 0.3*fermi;
500 e2 = 0.35*fermi;
501 }
502 else if (fParticle == theNeutron)
503 {
504 diffuse = 0.63*fermi;
505 // diffuse = 0.6*fermi;
506 G4double k0 = 1*GeV/hbarc;
507 diffuse *= k0/fWaveVector;
508 gamma = 0.3*fermi;
509 delta = 0.1*fermi*fermi;
510 e1 = 0.3*fermi;
511 e2 = 0.35*fermi;
512 }
513 else // as proton, if were not defined
514 {
515 diffuse = 0.63*fermi;
516 gamma = 0.3*fermi;
517 delta = 0.1*fermi*fermi;
518 e1 = 0.3*fermi;
519 e2 = 0.35*fermi;
520 }
521 G4double lambda = 15.; // 15 ok
522 // G4double kgamma = fWaveVector*gamma; // wavek*delta;
523 G4double kgamma = lambda*(1.-G4Exp(-fWaveVector*gamma/lambda)); // wavek*delta;
524
525 // G4cout<<"kgamma = "<<kgamma<<G4endl;
526
527 if(fAddCoulomb) // add Coulomb correction
528 {
529 G4double sinHalfTheta = std::sin(0.5*theta);
530 G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
531
532 kgamma += 0.5*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
533 // kgamma += 0.65*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
534 }
535
536 G4double kgamma2 = kgamma*kgamma;
537
538 // G4double dk2t = delta*fWaveVector*fWaveVector*theta; // delta*wavek*wavek*theta;
539 // G4cout<<"dk2t = "<<dk2t<<G4endl;
540 // G4double dk2t2 = dk2t*dk2t;
541 // G4double pikdt = pi*fWaveVector*diffuse*theta;// pi*wavek*diffuse*theta;
542
543 G4double pikdt = lambda*(1.-G4Exp(-pi*fWaveVector*diffuse*theta/lambda)); // wavek*delta;
544
545 // G4cout<<"pikdt = "<<pikdt<<G4endl;
546
547 damp = DampFactor(pikdt);
548 damp2 = damp*damp;
549
550 G4double mode2k2 = (e1*e1+e2*e2)*fWaveVector*fWaveVector;
551 G4double e2dk3t = -2.*e2*delta*fWaveVector*fWaveVector*fWaveVector*theta;
552
553 sigma = kgamma2;
554 // sigma += dk2t2;
555 sigma *= bzero2;
556 sigma += mode2k2*bone2;
557 sigma += e2dk3t*bzero*bone;
558
559 // sigma += kr2*(1 + 8.*fZommerfeld*fZommerfeld/kr2)*bonebyarg2; // correction at J1()/()
560 sigma += kr2*bonebyarg2; // correction at J1()/()
561
562 sigma *= damp2; // *rad*rad;
563
564 return sigma;
565}

Referenced by GetDiffuseElasticSumXsc().

◆ GetDiffElasticSumProbA()

G4double G4DiffuseElastic::GetDiffElasticSumProbA ( G4double alpha)

Definition at line 574 of file G4DiffuseElastic.cc.

575{
576 G4double theta;
577
578 theta = std::sqrt(alpha);
579
580 // theta = std::acos( 1 - alpha/2. );
581
582 G4double sigma, bzero, bzero2, bonebyarg, bonebyarg2, damp, damp2;
583 G4double delta, diffuse, gamma;
584 G4double e1, e2, bone, bone2;
585
586 // G4double wavek = momentum/hbarc; // wave vector
587 // G4double r0 = 1.08*fermi;
588 // G4double rad = r0*G4Pow::GetInstance()->A13(A);
589
590 G4double kr = fWaveVector*fNuclearRadius; // wavek*rad;
591 G4double kr2 = kr*kr;
592 G4double krt = kr*theta;
593
594 bzero = BesselJzero(krt);
595 bzero2 = bzero*bzero;
596 bone = BesselJone(krt);
597 bone2 = bone*bone;
598 bonebyarg = BesselOneByArg(krt);
599 bonebyarg2 = bonebyarg*bonebyarg;
600
601 if ( fParticle == theProton )
602 {
603 diffuse = 0.63*fermi;
604 // diffuse = 0.6*fermi;
605 gamma = 0.3*fermi;
606 delta = 0.1*fermi*fermi;
607 e1 = 0.3*fermi;
608 e2 = 0.35*fermi;
609 }
610 else if ( fParticle == theNeutron )
611 {
612 diffuse = 0.63*fermi;
613 // diffuse = 0.6*fermi;
614 // G4double k0 = 0.8*GeV/hbarc;
615 // diffuse *= k0/fWaveVector;
616 gamma = 0.3*fermi;
617 delta = 0.1*fermi*fermi;
618 e1 = 0.3*fermi;
619 e2 = 0.35*fermi;
620 }
621 else // as proton, if were not defined
622 {
623 diffuse = 0.63*fermi;
624 gamma = 0.3*fermi;
625 delta = 0.1*fermi*fermi;
626 e1 = 0.3*fermi;
627 e2 = 0.35*fermi;
628 }
629 G4double lambda = 15.; // 15 ok
630 // G4double kgamma = fWaveVector*gamma; // wavek*delta;
631 G4double kgamma = lambda*(1.-G4Exp(-fWaveVector*gamma/lambda)); // wavek*delta;
632
633 // G4cout<<"kgamma = "<<kgamma<<G4endl;
634
635 if( fAddCoulomb ) // add Coulomb correction
636 {
637 G4double sinHalfTheta = theta*0.5; // std::sin(0.5*theta);
638 G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
639
640 kgamma += 0.5*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
641 // kgamma += 0.65*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
642 }
643 G4double kgamma2 = kgamma*kgamma;
644
645 // G4double dk2t = delta*fWaveVector*fWaveVector*theta; // delta*wavek*wavek*theta;
646 // G4cout<<"dk2t = "<<dk2t<<G4endl;
647 // G4double dk2t2 = dk2t*dk2t;
648 // G4double pikdt = pi*fWaveVector*diffuse*theta;// pi*wavek*diffuse*theta;
649
650 G4double pikdt = lambda*(1. - G4Exp( -pi*fWaveVector*diffuse*theta/lambda ) ); // wavek*delta;
651
652 // G4cout<<"pikdt = "<<pikdt<<G4endl;
653
654 damp = DampFactor( pikdt );
655 damp2 = damp*damp;
656
657 G4double mode2k2 = ( e1*e1 + e2*e2 )*fWaveVector*fWaveVector;
658 G4double e2dk3t = -2.*e2*delta*fWaveVector*fWaveVector*fWaveVector*theta;
659
660 sigma = kgamma2;
661 // sigma += dk2t2;
662 sigma *= bzero2;
663 sigma += mode2k2*bone2;
664 sigma += e2dk3t*bzero*bone;
665
666 // sigma += kr2*(1 + 8.*fZommerfeld*fZommerfeld/kr2)*bonebyarg2; // correction at J1()/()
667 sigma += kr2*bonebyarg2; // correction at J1()/()
668
669 sigma *= damp2; // *rad*rad;
670
671 return sigma;
672}

Referenced by GetIntegrandFunction().

◆ GetDiffuseElasticSumXsc()

G4double G4DiffuseElastic::GetDiffuseElasticSumXsc ( const G4ParticleDefinition * particle,
G4double theta,
G4double momentum,
G4double A,
G4double Z )

Definition at line 244 of file G4DiffuseElastic.cc.

248{
249 fParticle = particle;
250 fWaveVector = momentum/hbarc;
251 fAtomicWeight = A;
252 fAtomicNumber = Z;
253 fNuclearRadius = CalculateNuclearRad(A);
254 fAddCoulomb = false;
255
256 G4double z = particle->GetPDGCharge();
257
258 G4double kRt = fWaveVector*fNuclearRadius*theta;
259 G4double kRtC = 1.9;
260
261 if( z && (kRt > kRtC) )
262 {
263 fAddCoulomb = true;
264 fBeta = CalculateParticleBeta( particle, momentum);
265 fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
266 fAm = CalculateAm( momentum, fZommerfeld, fAtomicNumber);
267 }
268 G4double sigma = fNuclearRadius*fNuclearRadius*GetDiffElasticSumProb(theta);
269
270 return sigma;
271}
G4double GetDiffElasticSumProb(G4double theta)
G4double CalculateNuclearRad(G4double A)

Referenced by GetInvElasticSumXsc().

◆ GetDiffuseElasticXsc()

G4double G4DiffuseElastic::GetDiffuseElasticXsc ( const G4ParticleDefinition * particle,
G4double theta,
G4double momentum,
G4double A )

Definition at line 173 of file G4DiffuseElastic.cc.

177{
178 fParticle = particle;
179 fWaveVector = momentum/hbarc;
180 fAtomicWeight = A;
181 fAddCoulomb = false;
182 fNuclearRadius = CalculateNuclearRad(A);
183
184 G4double sigma = fNuclearRadius*fNuclearRadius*GetDiffElasticProb(theta);
185
186 return sigma;
187}
G4double GetDiffElasticProb(G4double theta)

Referenced by GetInvElasticXsc().

◆ GetIntegrandFunction()

G4double G4DiffuseElastic::GetIntegrandFunction ( G4double theta)

Definition at line 680 of file G4DiffuseElastic.cc.

681{
682 G4double result;
683
684 result = GetDiffElasticSumProbA(alpha);
685
686 // result *= 2*pi*std::sin(theta);
687
688 return result;
689}
G4double GetDiffElasticSumProbA(G4double alpha)

Referenced by BuildAngleTable(), IntegralElasticProb(), SampleThetaCMS(), and TestAngleTable().

◆ GetInvCoulombElasticXsc()

G4double G4DiffuseElastic::GetInvCoulombElasticXsc ( const G4ParticleDefinition * particle,
G4double tMand,
G4double momentum,
G4double A,
G4double Z )

Definition at line 331 of file G4DiffuseElastic.cc.

335{
336 G4double m1 = particle->GetPDGMass();
337 G4LorentzVector lv1(0.,0.,plab,std::sqrt(plab*plab+m1*m1));
338
339 G4int iZ = static_cast<G4int>(Z+0.5);
340 G4int iA = static_cast<G4int>(A+0.5);
341 G4ParticleDefinition * theDef = 0;
342
343 if (iZ == 1 && iA == 1) theDef = theProton;
344 else if (iZ == 1 && iA == 2) theDef = theDeuteron;
345 else if (iZ == 1 && iA == 3) theDef = G4Triton::Triton();
346 else if (iZ == 2 && iA == 3) theDef = G4He3::He3();
347 else if (iZ == 2 && iA == 4) theDef = theAlpha;
348 else theDef = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(iZ,iA,0);
349
350 G4double tmass = theDef->GetPDGMass();
351
352 G4LorentzVector lv(0.0,0.0,0.0,tmass);
353 lv += lv1;
354
355 G4ThreeVector bst = lv.boostVector();
356 lv1.boost(-bst);
357
358 G4ThreeVector p1 = lv1.vect();
359 G4double ptot = p1.mag();
360 G4double ptot2 = ptot*ptot;
361 G4double cost = 1 - 0.5*std::fabs(tMand)/ptot2;
362
363 if( cost >= 1.0 ) cost = 1.0;
364 else if( cost <= -1.0) cost = -1.0;
365
366 G4double thetaCMS = std::acos(cost);
367
368 G4double sigma = GetCoulombElasticXsc( particle, thetaCMS, ptot, Z );
369
370 sigma *= pi/ptot2;
371
372 return sigma;
373}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double mag() const
G4double GetCoulombElasticXsc(const G4ParticleDefinition *particle, G4double theta, G4double momentum, G4double Z)
static G4He3 * He3()
Definition G4He3.cc:90
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Triton * Triton()
Definition G4Triton.cc:90
const G4double pi

◆ GetInvElasticSumXsc()

G4double G4DiffuseElastic::GetInvElasticSumXsc ( const G4ParticleDefinition * particle,
G4double tMand,
G4double momentum,
G4double A,
G4double Z )

Definition at line 279 of file G4DiffuseElastic.cc.

283{
284 G4double m1 = particle->GetPDGMass();
285
286 G4LorentzVector lv1(0.,0.,plab,std::sqrt(plab*plab+m1*m1));
287
288 G4int iZ = static_cast<G4int>(Z+0.5);
289 G4int iA = static_cast<G4int>(A+0.5);
290
291 G4ParticleDefinition* theDef = 0;
292
293 if (iZ == 1 && iA == 1) theDef = theProton;
294 else if (iZ == 1 && iA == 2) theDef = theDeuteron;
295 else if (iZ == 1 && iA == 3) theDef = G4Triton::Triton();
296 else if (iZ == 2 && iA == 3) theDef = G4He3::He3();
297 else if (iZ == 2 && iA == 4) theDef = theAlpha;
298 else theDef = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(iZ,iA,0);
299
300 G4double tmass = theDef->GetPDGMass();
301
302 G4LorentzVector lv(0.0,0.0,0.0,tmass);
303 lv += lv1;
304
305 G4ThreeVector bst = lv.boostVector();
306 lv1.boost(-bst);
307
308 G4ThreeVector p1 = lv1.vect();
309 G4double ptot = p1.mag();
310 G4double ptot2 = ptot*ptot;
311 G4double cost = 1 - 0.5*std::fabs(tMand)/ptot2;
312
313 if( cost >= 1.0 ) cost = 1.0;
314 else if( cost <= -1.0) cost = -1.0;
315
316 G4double thetaCMS = std::acos(cost);
317
318 G4double sigma = GetDiffuseElasticSumXsc( particle, thetaCMS, ptot, A, Z );
319
320 sigma *= pi/ptot2;
321
322 return sigma;
323}
G4double GetDiffuseElasticSumXsc(const G4ParticleDefinition *particle, G4double theta, G4double momentum, G4double A, G4double Z)

◆ GetInvElasticXsc()

G4double G4DiffuseElastic::GetInvElasticXsc ( const G4ParticleDefinition * particle,
G4double theta,
G4double momentum,
G4double A,
G4double Z )

Definition at line 194 of file G4DiffuseElastic.cc.

198{
199 G4double m1 = particle->GetPDGMass();
200 G4LorentzVector lv1(0.,0.,plab,std::sqrt(plab*plab+m1*m1));
201
202 G4int iZ = static_cast<G4int>(Z+0.5);
203 G4int iA = static_cast<G4int>(A+0.5);
204 G4ParticleDefinition * theDef = 0;
205
206 if (iZ == 1 && iA == 1) theDef = theProton;
207 else if (iZ == 1 && iA == 2) theDef = theDeuteron;
208 else if (iZ == 1 && iA == 3) theDef = G4Triton::Triton();
209 else if (iZ == 2 && iA == 3) theDef = G4He3::He3();
210 else if (iZ == 2 && iA == 4) theDef = theAlpha;
211 else theDef = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(iZ,iA,0);
212
213 G4double tmass = theDef->GetPDGMass();
214
215 G4LorentzVector lv(0.0,0.0,0.0,tmass);
216 lv += lv1;
217
218 G4ThreeVector bst = lv.boostVector();
219 lv1.boost(-bst);
220
221 G4ThreeVector p1 = lv1.vect();
222 G4double ptot = p1.mag();
223 G4double ptot2 = ptot*ptot;
224 G4double cost = 1 - 0.5*std::fabs(tMand)/ptot2;
225
226 if( cost >= 1.0 ) cost = 1.0;
227 else if( cost <= -1.0) cost = -1.0;
228
229 G4double thetaCMS = std::acos(cost);
230
231 G4double sigma = GetDiffuseElasticXsc( particle, thetaCMS, ptot, A);
232
233 sigma *= pi/ptot2;
234
235 return sigma;
236}
G4double GetDiffuseElasticXsc(const G4ParticleDefinition *particle, G4double theta, G4double momentum, G4double A)

◆ GetNuclearRadius()

G4double G4DiffuseElastic::GetNuclearRadius ( )
inline

Definition at line 191 of file G4DiffuseElastic.hh.

191{return fNuclearRadius;};

◆ GetScatteringAngle()

G4double G4DiffuseElastic::GetScatteringAngle ( G4int iMomentum,
G4int iAngle,
G4double position )

Definition at line 1094 of file G4DiffuseElastic.cc.

1095{
1096 G4double x1, x2, y1, y2, randAngle;
1097
1098 if( iAngle == 0 )
1099 {
1100 randAngle = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle);
1101 // iAngle++;
1102 }
1103 else
1104 {
1105 if ( iAngle >= G4int((*fAngleTable)(iMomentum)->GetVectorLength()) )
1106 {
1107 iAngle = G4int((*fAngleTable)(iMomentum)->GetVectorLength() - 1);
1108 }
1109 y1 = (*(*fAngleTable)(iMomentum))(iAngle-1);
1110 y2 = (*(*fAngleTable)(iMomentum))(iAngle);
1111
1112 x1 = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle-1);
1113 x2 = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle);
1114
1115 if ( x1 == x2 ) randAngle = x2;
1116 else
1117 {
1118 if ( y1 == y2 ) randAngle = x1 + ( x2 - x1 )*G4UniformRand();
1119 else
1120 {
1121 randAngle = x1 + ( position - y1 )*( x2 - x1 )/( y2 - y1 );
1122 }
1123 }
1124 }
1125 return randAngle;
1126}
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by SampleTableThetaCMS().

◆ Initialise()

void G4DiffuseElastic::Initialise ( )

Definition at line 139 of file G4DiffuseElastic.cc.

140{
141
142 // fEnergyVector = new G4PhysicsLogVector( theMinEnergy, theMaxEnergy, fEnergyBin );
143
144 const G4ElementTable* theElementTable = G4Element::GetElementTable();
145
146 std::size_t jEl, numOfEl = G4Element::GetNumberOfElements();
147
148 for( jEl = 0; jEl < numOfEl; ++jEl) // application element loop
149 {
150 fAtomicNumber = (*theElementTable)[jEl]->GetZ(); // atomic number
151 fAtomicWeight = G4NistManager::Instance()->GetAtomicMassAmu( static_cast< G4int >( fAtomicNumber ) );
152 fNuclearRadius = CalculateNuclearRad(fAtomicWeight);
153
154 if( verboseLevel > 0 )
155 {
156 G4cout<<"G4DiffuseElastic::Initialise() the element: "
157 <<(*theElementTable)[jEl]->GetName()<<G4endl;
158 }
159 fElementNumberVector.push_back(fAtomicNumber);
160 fElementNameVector.push_back((*theElementTable)[jEl]->GetName());
161
163 fAngleBank.push_back(fAngleTable);
164 }
165 return;
166}
std::vector< G4Element * > G4ElementTable
static std::size_t GetNumberOfElements()
Definition G4Element.cc:408
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const

◆ InitialiseOnFly()

void G4DiffuseElastic::InitialiseOnFly ( G4double Z,
G4double A )

Definition at line 975 of file G4DiffuseElastic.cc.

976{
977 fAtomicNumber = Z; // atomic number
978 fAtomicWeight = G4NistManager::Instance()->GetAtomicMassAmu( static_cast< G4int >( Z ) );
979
980 fNuclearRadius = CalculateNuclearRad(fAtomicWeight);
981
982 if( verboseLevel > 0 )
983 {
984 G4cout<<"G4DiffuseElastic::InitialiseOnFly() the element with Z = "
985 <<Z<<"; and A = "<<A<<G4endl;
986 }
987 fElementNumberVector.push_back(fAtomicNumber);
988
990
991 fAngleBank.push_back(fAngleTable);
992
993 return;
994}

Referenced by SampleTableThetaCMS().

◆ IntegralElasticProb()

G4double G4DiffuseElastic::IntegralElasticProb ( const G4ParticleDefinition * particle,
G4double theta,
G4double momentum,
G4double A )

Definition at line 696 of file G4DiffuseElastic.cc.

700{
701 G4double result;
702 fParticle = particle;
703 fWaveVector = momentum/hbarc;
704 fAtomicWeight = A;
705
706 fNuclearRadius = CalculateNuclearRad(A);
707
708
709 G4Integrator<G4DiffuseElastic,G4double(G4DiffuseElastic::*)(G4double)> integral;
710
711 // result = integral.Legendre10(this,&G4DiffuseElastic::GetIntegrandFunction, 0., theta );
712 result = integral.Legendre96(this,&G4DiffuseElastic::GetIntegrandFunction, 0., theta );
713
714 return result;
715}

◆ IsApplicable()

G4bool G4DiffuseElastic::IsApplicable ( const G4HadProjectile & projectile,
G4Nucleus & nucleus )
inlinevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 232 of file G4DiffuseElastic.hh.

234{
235 if( ( projectile.GetDefinition() == G4Proton::Proton() ||
236 projectile.GetDefinition() == G4Neutron::Neutron() ||
237 projectile.GetDefinition() == G4PionPlus::PionPlus() ||
238 projectile.GetDefinition() == G4PionMinus::PionMinus() ||
239 projectile.GetDefinition() == G4KaonPlus::KaonPlus() ||
240 projectile.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
241
242 nucleus.GetZ_asInt() >= 2 ) return true;
243 else return false;
244}
const G4ParticleDefinition * GetDefinition() const
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()

◆ NeutronTuniform()

G4double G4DiffuseElastic::NeutronTuniform ( G4int Z)

Definition at line 826 of file G4DiffuseElastic.cc.

827{
828 G4double elZ = G4double(Z);
829 elZ -= 1.;
830 // G4double Tkin = 20.*G4Exp(-elZ/10.) + 1.;
831 G4double Tkin = 12.*G4Exp(-elZ/10.) + 1.;
832 return Tkin;
833}

Referenced by SampleInvariantT().

◆ SampleInvariantT()

G4double G4DiffuseElastic::SampleInvariantT ( const G4ParticleDefinition * p,
G4double plab,
G4int Z,
G4int A )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 788 of file G4DiffuseElastic.cc.

790{
791 fParticle = aParticle;
792 G4double m1 = fParticle->GetPDGMass(), t;
793 G4double totElab = std::sqrt(m1*m1+p*p);
795 G4LorentzVector lv1(p,0.0,0.0,totElab);
796 G4LorentzVector lv(0.0,0.0,0.0,mass2);
797 lv += lv1;
798
799 G4ThreeVector bst = lv.boostVector();
800 lv1.boost(-bst);
801
802 G4ThreeVector p1 = lv1.vect();
803 G4double momentumCMS = p1.mag();
804
805 if( aParticle == theNeutron)
806 {
807 G4double Tmax = NeutronTuniform( Z );
808 G4double pCMS2 = momentumCMS*momentumCMS;
809 G4double Tkin = std::sqrt(pCMS2+m1*m1)-m1;
810
811 if( Tkin <= Tmax )
812 {
813 t = 4.*pCMS2*G4UniformRand();
814 // G4cout<<Tkin<<", "<<Tmax<<", "<<std::sqrt(t)<<"; ";
815 return t;
816 }
817 }
818
819 t = SampleTableT( aParticle, momentumCMS, G4double(Z), G4double(A) ); // sample theta2 in cms
820
821 return t;
822}
G4double SampleTableT(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)
G4double NeutronTuniform(G4int Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ SampleT()

G4double G4DiffuseElastic::SampleT ( const G4ParticleDefinition * aParticle,
G4double p,
G4double A )

Definition at line 721 of file G4DiffuseElastic.cc.

722{
723 G4double theta = SampleThetaCMS( aParticle, p, A); // sample theta in cms
724 G4double t = 2*p*p*( 1 - std::cos(theta) ); // -t !!!
725 return t;
726}
G4double SampleThetaCMS(const G4ParticleDefinition *aParticle, G4double p, G4double A)

Referenced by SampleThetaLab().

◆ SampleTableT()

G4double G4DiffuseElastic::SampleTableT ( const G4ParticleDefinition * aParticle,
G4double p,
G4double Z,
G4double A )

Definition at line 840 of file G4DiffuseElastic.cc.

842{
843 G4double alpha = SampleTableThetaCMS( aParticle, p, Z, A); // sample theta2 in cms
844 G4double t = 2*p*p*( 1 - std::cos(std::sqrt(alpha)) ); // -t !!!
845 // G4double t = p*p*alpha; // -t !!!
846 return t;
847}
G4double SampleTableThetaCMS(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)

Referenced by SampleInvariantT().

◆ SampleTableThetaCMS()

G4double G4DiffuseElastic::SampleTableThetaCMS ( const G4ParticleDefinition * aParticle,
G4double p,
G4double Z,
G4double A )

Definition at line 855 of file G4DiffuseElastic.cc.

857{
858 std::size_t iElement;
859 G4int iMomentum, iAngle;
860 G4double randAngle, position, theta1, theta2, E1, E2, W1, W2, W;
861 G4double m1 = particle->GetPDGMass();
862
863 for(iElement = 0; iElement < fElementNumberVector.size(); iElement++)
864 {
865 if( std::fabs(Z - fElementNumberVector[iElement]) < 0.5) break;
866 }
867 if ( iElement == fElementNumberVector.size() )
868 {
869 InitialiseOnFly(Z,A); // table preparation, if needed
870
871 // iElement--;
872
873 // G4cout << "G4DiffuseElastic: Element with atomic number " << Z
874 // << " is not found, return zero angle" << G4endl;
875 // return 0.; // no table for this element
876 }
877 // G4cout<<"iElement = "<<iElement<<G4endl;
878
879 fAngleTable = fAngleBank[iElement];
880
881 G4double kinE = std::sqrt(momentum*momentum + m1*m1) - m1;
882
883 for( iMomentum = 0; iMomentum < fEnergyBin; iMomentum++)
884 {
885 if( kinE < fEnergyVector->GetLowEdgeEnergy(iMomentum) ) break;
886 }
887 if ( iMomentum >= fEnergyBin ) iMomentum = fEnergyBin-1; // kinE is more then theMaxEnergy
888 if ( iMomentum < 0 ) iMomentum = 0; // against negative index, kinE < theMinEnergy
889
890 // G4cout<<"iMomentum = "<<iMomentum<<G4endl;
891
892 if (iMomentum == fEnergyBin -1 || iMomentum == 0 ) // the table edges
893 {
894 position = (*(*fAngleTable)(iMomentum))(fAngleBin-2)*G4UniformRand();
895
896 // G4cout<<"position = "<<position<<G4endl;
897
898 for(iAngle = 0; iAngle < fAngleBin-1; iAngle++)
899 {
900 if( position > (*(*fAngleTable)(iMomentum))(iAngle) ) break;
901 }
902 if (iAngle >= fAngleBin-1) iAngle = fAngleBin-2;
903
904 // G4cout<<"iAngle = "<<iAngle<<G4endl;
905
906 randAngle = GetScatteringAngle(iMomentum, iAngle, position);
907
908 // G4cout<<"randAngle = "<<randAngle<<G4endl;
909 }
910 else // kinE inside between energy table edges
911 {
912 // position = (*(*fAngleTable)(iMomentum))(fAngleBin-2)*G4UniformRand();
913 position = (*(*fAngleTable)(iMomentum))(0)*G4UniformRand();
914
915 // G4cout<<"position = "<<position<<G4endl;
916
917 for(iAngle = 0; iAngle < fAngleBin-1; iAngle++)
918 {
919 // if( position < (*(*fAngleTable)(iMomentum))(iAngle) ) break;
920 if( position > (*(*fAngleTable)(iMomentum))(iAngle) ) break;
921 }
922 if (iAngle >= fAngleBin-1) iAngle = fAngleBin-2;
923
924 // G4cout<<"iAngle = "<<iAngle<<G4endl;
925
926 theta2 = GetScatteringAngle(iMomentum, iAngle, position);
927
928 // G4cout<<"theta2 = "<<theta2<<G4endl;
929 E2 = fEnergyVector->GetLowEdgeEnergy(iMomentum);
930
931 // G4cout<<"E2 = "<<E2<<G4endl;
932
933 iMomentum--;
934
935 // position = (*(*fAngleTable)(iMomentum))(fAngleBin-2)*G4UniformRand();
936
937 // G4cout<<"position = "<<position<<G4endl;
938
939 for(iAngle = 0; iAngle < fAngleBin-1; iAngle++)
940 {
941 // if( position < (*(*fAngleTable)(iMomentum))(iAngle) ) break;
942 if( position > (*(*fAngleTable)(iMomentum))(iAngle) ) break;
943 }
944 if (iAngle >= fAngleBin-1) iAngle = fAngleBin-2;
945
946 theta1 = GetScatteringAngle(iMomentum, iAngle, position);
947
948 // G4cout<<"theta1 = "<<theta1<<G4endl;
949
950 E1 = fEnergyVector->GetLowEdgeEnergy(iMomentum);
951
952 // G4cout<<"E1 = "<<E1<<G4endl;
953
954 W = 1.0/(E2 - E1);
955 W1 = (E2 - kinE)*W;
956 W2 = (kinE - E1)*W;
957
958 randAngle = W1*theta1 + W2*theta2;
959
960 // randAngle = theta2;
961 // G4cout<<"randAngle = "<<randAngle<<G4endl;
962 }
963 // G4double angle = randAngle;
964 // if (randAngle > 0.) randAngle /= 2*pi*std::sin(angle);
965
966 if(randAngle < 0.) randAngle = 0.;
967
968 return randAngle;
969}
void InitialiseOnFly(G4double Z, G4double A)
G4double GetScatteringAngle(G4int iMomentum, G4int iAngle, G4double position)
#define W
Definition crc32.c:85

Referenced by SampleTableT().

◆ SampleThetaCMS()

G4double G4DiffuseElastic::SampleThetaCMS ( const G4ParticleDefinition * aParticle,
G4double p,
G4double A )

Definition at line 734 of file G4DiffuseElastic.cc.

736{
737 G4int i, iMax = 100;
738 G4double norm, theta1, theta2, thetaMax;
739 G4double result = 0., sum = 0.;
740
741 fParticle = particle;
742 fWaveVector = momentum/hbarc;
743 fAtomicWeight = A;
744
745 fNuclearRadius = CalculateNuclearRad(A);
746
747 thetaMax = 10.174/fWaveVector/fNuclearRadius;
748
749 if (thetaMax > pi) thetaMax = pi;
750
751 G4Integrator<G4DiffuseElastic,G4double(G4DiffuseElastic::*)(G4double)> integral;
752
753 // result = integral.Legendre10(this,&G4DiffuseElastic::GetIntegrandFunction, 0., theta );
754 norm = integral.Legendre96(this,&G4DiffuseElastic::GetIntegrandFunction, 0., thetaMax );
755
756 norm *= G4UniformRand();
757
758 for(i = 1; i <= iMax; i++)
759 {
760 theta1 = (i-1)*thetaMax/iMax;
761 theta2 = i*thetaMax/iMax;
762 sum += integral.Legendre10(this,&G4DiffuseElastic::GetIntegrandFunction, theta1, theta2);
763
764 if ( sum >= norm )
765 {
766 result = 0.5*(theta1 + theta2);
767 break;
768 }
769 }
770 if (i > iMax ) result = 0.5*(theta1 + theta2);
771
772 G4double sigma = pi*thetaMax/iMax;
773
774 result += G4RandGauss::shoot(0.,sigma);
775
776 if(result < 0.) result = 0.;
777 if(result > thetaMax) result = thetaMax;
778
779 return result;
780}

Referenced by SampleT().

◆ SampleThetaLab()

G4double G4DiffuseElastic::SampleThetaLab ( const G4HadProjectile * aParticle,
G4double tmass,
G4double A )

Definition at line 1137 of file G4DiffuseElastic.cc.

1139{
1140 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1141 G4double m1 = theParticle->GetPDGMass();
1142 G4double plab = aParticle->GetTotalMomentum();
1143 G4LorentzVector lv1 = aParticle->Get4Momentum();
1144 G4LorentzVector lv(0.0,0.0,0.0,tmass);
1145 lv += lv1;
1146
1147 G4ThreeVector bst = lv.boostVector();
1148 lv1.boost(-bst);
1149
1150 G4ThreeVector p1 = lv1.vect();
1151 G4double ptot = p1.mag();
1152 G4double tmax = 4.0*ptot*ptot;
1153 G4double t = 0.0;
1154
1155
1156 //
1157 // Sample t
1158 //
1159
1160 t = SampleT( theParticle, ptot, A);
1161 if (t <= 0.0) { return 0.0; }
1162
1163 // NaN finder substituted
1164 if (t > tmax)
1165 {
1166 if (verboseLevel > 0)
1167 {
1168 G4cout << "G4DiffuseElastic:WARNING: A = " << A
1169 << " mom(GeV)= " << plab/GeV
1170 << " S-wave will be sampled"
1171 << G4endl;
1172 }
1173 t = G4UniformRand()*tmax;
1174 }
1175 if(verboseLevel>1)
1176 {
1177 G4cout <<" t= " << t << " tmax= " << tmax
1178 << " ptot= " << ptot << G4endl;
1179 }
1180 // Sampling of angles in CM system
1181
1182 G4double phi = G4UniformRand()*twopi;
1183 G4double cost = 1. - 2.0*t/tmax;
1184 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
1185
1186 if (verboseLevel>1)
1187 {
1188 G4cout << "cos(t)=" << cost << " std::sin(t)=" << sint << G4endl;
1189 }
1190 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
1191 v1 *= ptot;
1192 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1));
1193
1194 nlv1.boost(bst);
1195
1196 G4ThreeVector np1 = nlv1.vect();
1197
1198 // G4double theta = std::acos( np1.z()/np1.mag() ); // degree;
1199
1200 G4double theta = np1.theta();
1201
1202 return theta;
1203}
double theta() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
G4double SampleT(const G4ParticleDefinition *aParticle, G4double p, G4double A)
G4double GetTotalMomentum() const
const G4LorentzVector & Get4Momentum() const

◆ SetHEModelLowLimit()

void G4DiffuseElastic::SetHEModelLowLimit ( G4double value)
inline

Definition at line 256 of file G4DiffuseElastic.hh.

257{
258 lowEnergyLimitHE = value;
259}

◆ SetLowestEnergyLimit()

void G4DiffuseElastic::SetLowestEnergyLimit ( G4double value)
inline

Definition at line 266 of file G4DiffuseElastic.hh.

267{
268 lowestEnergyLimit = value;
269}

◆ SetPlabLowLimit()

void G4DiffuseElastic::SetPlabLowLimit ( G4double value)
inline

Definition at line 251 of file G4DiffuseElastic.hh.

252{
253 plabLowLimit = value;
254}

◆ SetQModelLowLimit()

void G4DiffuseElastic::SetQModelLowLimit ( G4double value)
inline

Definition at line 261 of file G4DiffuseElastic.hh.

262{
263 lowEnergyLimitQ = value;
264}

◆ SetRecoilKinEnergyLimit()

void G4DiffuseElastic::SetRecoilKinEnergyLimit ( G4double value)
inline

Definition at line 246 of file G4DiffuseElastic.hh.

247{
248 lowEnergyRecoilLimit = value;
249}

◆ TestAngleTable()

void G4DiffuseElastic::TestAngleTable ( const G4ParticleDefinition * theParticle,
G4double partMom,
G4double Z,
G4double A )

Definition at line 1331 of file G4DiffuseElastic.cc.

1333{
1334 fAtomicNumber = Z; // atomic number
1335 fAtomicWeight = A; // number of nucleons
1336 fNuclearRadius = CalculateNuclearRad(fAtomicWeight);
1337
1338
1339
1340 G4cout<<"G4DiffuseElastic::TestAngleTable() init the element with Z = "
1341 <<Z<<"; and A = "<<A<<G4endl;
1342
1343 fElementNumberVector.push_back(fAtomicNumber);
1344
1345
1346
1347
1348 G4int i=0, j;
1349 G4double a = 0., z = theParticle->GetPDGCharge(), m1 = fParticle->GetPDGMass();
1350 G4double alpha1=0., alpha2=0., alphaMax=0., alphaCoulomb=0.;
1351 G4double deltaL10 = 0., deltaL96 = 0., deltaAG = 0.;
1352 G4double sumL10 = 0.,sumL96 = 0.,sumAG = 0.;
1353 G4double epsilon = 0.001;
1354
1355 G4Integrator<G4DiffuseElastic,G4double(G4DiffuseElastic::*)(G4double)> integral;
1356
1357 fAngleTable = new G4PhysicsTable(fEnergyBin);
1358
1359 fWaveVector = partMom/hbarc;
1360
1361 G4double kR = fWaveVector*fNuclearRadius;
1362 G4double kR2 = kR*kR;
1363 G4double kRmax = 10.6; // 10.6, 18, 10.174; ~ 3 maxima of J1 or 15., 25.
1364 G4double kRcoul = 1.2; // 1.4, 2.5; // on the first slope of J1
1365
1366 alphaMax = kRmax*kRmax/kR2;
1367
1368 if (alphaMax > 4.) alphaMax = 4.; // vmg05-02-09: was pi2
1369
1370 alphaCoulomb = kRcoul*kRcoul/kR2;
1371
1372 if( z )
1373 {
1374 a = partMom/m1; // beta*gamma for m1
1375 fBeta = a/std::sqrt(1+a*a);
1376 fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
1377 fAm = CalculateAm( partMom, fZommerfeld, fAtomicNumber);
1378 }
1379 G4PhysicsFreeVector* angleVector = new G4PhysicsFreeVector(fAngleBin-1);
1380
1381 // G4PhysicsLogVector* angleBins = new G4PhysicsLogVector( 0.001*alphaMax, alphaMax, fAngleBin );
1382
1383
1384 fAddCoulomb = false;
1385
1386 for(j = 1; j < fAngleBin; j++)
1387 {
1388 // alpha1 = angleBins->GetLowEdgeEnergy(j-1);
1389 // alpha2 = angleBins->GetLowEdgeEnergy(j);
1390
1391 alpha1 = alphaMax*(j-1)/fAngleBin;
1392 alpha2 = alphaMax*( j )/fAngleBin;
1393
1394 if( ( alpha2 > alphaCoulomb ) && z ) fAddCoulomb = true;
1395
1396 deltaL10 = integral.Legendre10(this, &G4DiffuseElastic::GetIntegrandFunction, alpha1, alpha2);
1397 deltaL96 = integral.Legendre96(this, &G4DiffuseElastic::GetIntegrandFunction, alpha1, alpha2);
1398 deltaAG = integral.AdaptiveGauss(this, &G4DiffuseElastic::GetIntegrandFunction,
1399 alpha1, alpha2,epsilon);
1400
1401 // G4cout<<alpha1<<"\t"<<std::sqrt(alpha1)/degree<<"\t"
1402 // <<deltaL10<<"\t"<<deltaL96<<"\t"<<deltaAG<<G4endl;
1403
1404 sumL10 += deltaL10;
1405 sumL96 += deltaL96;
1406 sumAG += deltaAG;
1407
1408 G4cout<<alpha1<<"\t"<<std::sqrt(alpha1)/degree<<"\t"
1409 <<sumL10<<"\t"<<sumL96<<"\t"<<sumAG<<G4endl;
1410
1411 angleVector->PutValue( j-1 , alpha1, sumL10 ); // alpha2
1412 }
1413 fAngleTable->insertAt(i,angleVector);
1414 fAngleBank.push_back(fAngleTable);
1415
1416 /*
1417 // Integral over all angle range - Bad accuracy !!!
1418
1419 sumL10 = integral.Legendre10(this, &G4DiffuseElastic::GetIntegrandFunction, 0., alpha2);
1420 sumL96 = integral.Legendre96(this, &G4DiffuseElastic::GetIntegrandFunction, 0., alpha2);
1421 sumAG = integral.AdaptiveGauss(this, &G4DiffuseElastic::GetIntegrandFunction,
1422 0., alpha2,epsilon);
1423 G4cout<<G4endl;
1424 G4cout<<alpha2<<"\t"<<std::sqrt(alpha2)/degree<<"\t"
1425 <<sumL10<<"\t"<<sumL96<<"\t"<<sumAG<<G4endl;
1426 */
1427 return;
1428}
G4double epsilon(G4double density, G4double temperature)

◆ ThetaCMStoThetaLab()

G4double G4DiffuseElastic::ThetaCMStoThetaLab ( const G4DynamicParticle * aParticle,
G4double tmass,
G4double thetaCMS )

Definition at line 1212 of file G4DiffuseElastic.cc.

1214{
1215 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1216 G4double m1 = theParticle->GetPDGMass();
1217 // G4double plab = aParticle->GetTotalMomentum();
1218 G4LorentzVector lv1 = aParticle->Get4Momentum();
1219 G4LorentzVector lv(0.0,0.0,0.0,tmass);
1220
1221 lv += lv1;
1222
1223 G4ThreeVector bst = lv.boostVector();
1224
1225 lv1.boost(-bst);
1226
1227 G4ThreeVector p1 = lv1.vect();
1228 G4double ptot = p1.mag();
1229
1230 G4double phi = G4UniformRand()*twopi;
1231 G4double cost = std::cos(thetaCMS);
1232 G4double sint;
1233
1234 if( cost >= 1.0 )
1235 {
1236 cost = 1.0;
1237 sint = 0.0;
1238 }
1239 else if( cost <= -1.0)
1240 {
1241 cost = -1.0;
1242 sint = 0.0;
1243 }
1244 else
1245 {
1246 sint = std::sqrt((1.0-cost)*(1.0+cost));
1247 }
1248 if (verboseLevel>1)
1249 {
1250 G4cout << "cos(tcms)=" << cost << " std::sin(tcms)=" << sint << G4endl;
1251 }
1252 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
1253 v1 *= ptot;
1254 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1));
1255
1256 nlv1.boost(bst);
1257
1258 G4ThreeVector np1 = nlv1.vect();
1259
1260
1261 G4double thetaLab = np1.theta();
1262
1263 return thetaLab;
1264}
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const

◆ ThetaLabToThetaCMS()

G4double G4DiffuseElastic::ThetaLabToThetaCMS ( const G4DynamicParticle * aParticle,
G4double tmass,
G4double thetaLab )

Definition at line 1272 of file G4DiffuseElastic.cc.

1274{
1275 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1276 G4double m1 = theParticle->GetPDGMass();
1277 G4double plab = aParticle->GetTotalMomentum();
1278 G4LorentzVector lv1 = aParticle->Get4Momentum();
1279 G4LorentzVector lv(0.0,0.0,0.0,tmass);
1280
1281 lv += lv1;
1282
1283 G4ThreeVector bst = lv.boostVector();
1284
1285 // lv1.boost(-bst);
1286
1287 // G4ThreeVector p1 = lv1.vect();
1288 // G4double ptot = p1.mag();
1289
1290 G4double phi = G4UniformRand()*twopi;
1291 G4double cost = std::cos(thetaLab);
1292 G4double sint;
1293
1294 if( cost >= 1.0 )
1295 {
1296 cost = 1.0;
1297 sint = 0.0;
1298 }
1299 else if( cost <= -1.0)
1300 {
1301 cost = -1.0;
1302 sint = 0.0;
1303 }
1304 else
1305 {
1306 sint = std::sqrt((1.0-cost)*(1.0+cost));
1307 }
1308 if (verboseLevel>1)
1309 {
1310 G4cout << "cos(tlab)=" << cost << " std::sin(tlab)=" << sint << G4endl;
1311 }
1312 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
1313 v1 *= plab;
1314 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(plab*plab + m1*m1));
1315
1316 nlv1.boost(-bst);
1317
1318 G4ThreeVector np1 = nlv1.vect();
1319
1320
1321 G4double thetaCMS = np1.theta();
1322
1323 return thetaCMS;
1324}
G4double GetTotalMomentum() const

The documentation for this class was generated from the following files: