Geant4 10.7.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 ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 ModelDescription (std::ostream &outFile) const
 
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
 
- 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:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Proton * Proton()
Definition: G4Proton.cc:92

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 1000 of file G4DiffuseElastic.cc.

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

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:120

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}
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:57
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:179
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}
double A(double temperature)
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}
double mag() const
G4double GetCoulombElasticXsc(const G4ParticleDefinition *particle, G4double theta, G4double momentum, G4double Z)
static G4He3 * He3()
Definition: G4He3.cc:93
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Triton * Triton()
Definition: G4Triton.cc:94
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 1093 of file G4DiffuseElastic.cc.

1094{
1095 G4double x1, x2, y1, y2, randAngle;
1096
1097 if( iAngle == 0 )
1098 {
1099 randAngle = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle);
1100 // iAngle++;
1101 }
1102 else
1103 {
1104 if ( iAngle >= G4int((*fAngleTable)(iMomentum)->GetVectorLength()) )
1105 {
1106 iAngle = (*fAngleTable)(iMomentum)->GetVectorLength() - 1;
1107 }
1108 y1 = (*(*fAngleTable)(iMomentum))(iAngle-1);
1109 y2 = (*(*fAngleTable)(iMomentum))(iAngle);
1110
1111 x1 = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle-1);
1112 x2 = (*fAngleTable)(iMomentum)->GetLowEdgeEnergy(iAngle);
1113
1114 if ( x1 == x2 ) randAngle = x2;
1115 else
1116 {
1117 if ( y1 == y2 ) randAngle = x1 + ( x2 - x1 )*G4UniformRand();
1118 else
1119 {
1120 randAngle = x1 + ( position - y1 )*( x2 - x1 )/( y2 - y1 );
1121 }
1122 }
1123 }
1124 return randAngle;
1125}
#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 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 G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const

◆ InitialiseOnFly()

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

Definition at line 974 of file G4DiffuseElastic.cc.

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

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
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()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115

◆ NeutronTuniform()

G4double G4DiffuseElastic::NeutronTuniform ( G4int  Z)

Definition at line 825 of file G4DiffuseElastic.cc.

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

Referenced by SampleInvariantT().

◆ SampleInvariantT()

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

Reimplemented from G4HadronicInteraction.

Definition at line 787 of file G4DiffuseElastic.cc.

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

841{
842 G4double alpha = SampleTableThetaCMS( aParticle, p, Z, A); // sample theta2 in cms
843 G4double t = 2*p*p*( 1 - std::cos(std::sqrt(alpha)) ); // -t !!!
844 // G4double t = p*p*alpha; // -t !!!
845 return t;
846}
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 854 of file G4DiffuseElastic.cc.

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

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

Referenced by SampleT().

◆ SampleThetaLab()

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

Definition at line 1136 of file G4DiffuseElastic.cc.

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

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

◆ ThetaCMStoThetaLab()

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

Definition at line 1224 of file G4DiffuseElastic.cc.

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

◆ ThetaLabToThetaCMS()

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

Definition at line 1284 of file G4DiffuseElastic.cc.

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

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