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

#include <G4DiffuseElasticV2.hh>

+ Inheritance diagram for G4DiffuseElasticV2:

Public Member Functions

 G4DiffuseElasticV2 ()
 
virtual ~G4DiffuseElasticV2 ()
 
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 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, unsigned long iAngle, G4double position)
 
G4double SampleThetaLab (const G4HadProjectile *aParticle, G4double tmass, G4double A)
 
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)
 
G4double BesselJzero (G4double z)
 
G4double BesselJone (G4double z)
 
G4double DampFactor (G4double z)
 
G4double BesselOneByArg (G4double z)
 
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 61 of file G4DiffuseElasticV2.hh.

Constructor & Destructor Documentation

◆ G4DiffuseElasticV2()

G4DiffuseElasticV2::G4DiffuseElasticV2 ( )

Definition at line 76 of file G4DiffuseElasticV2.cc.

77 : G4HadronElastic("DiffuseElasticV2"), fParticle(0)
78{
79 SetMinEnergy( 0.01*MeV );
81
82 verboseLevel = 0;
83 lowEnergyRecoilLimit = 100.*keV;
84 lowEnergyLimitQ = 0.0*GeV;
85 lowEnergyLimitHE = 0.0*GeV;
86 lowestEnergyLimit = 0.0*keV;
87 plabLowLimit = 20.0*MeV;
88
89 theProton = G4Proton::Proton();
90 theNeutron = G4Neutron::Neutron();
91
92 fEnergyBin = 300; // Increased from the original 200 to have no wider log-energy-bins up to 10 PeV
93 fAngleBin = 200;
94
95 fEnergyVector = new G4PhysicsLogVector( theMinEnergy, theMaxEnergy, fEnergyBin );
96
97 fEnergyAngleVector = 0;
98 fEnergySumVector = 0;
99
100 fParticle = 0;
101 fWaveVector = 0.;
102 fAtomicWeight = 0.;
103 fAtomicNumber = 0.;
104 fNuclearRadius = 0.;
105 fBeta = 0.;
106 fZommerfeld = 0.;
107 fAm = 0.;
108 fAddCoulomb = false;
109}
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 G4Proton * Proton()
Definition G4Proton.cc:90

Referenced by BuildAngleTable().

◆ ~G4DiffuseElasticV2()

G4DiffuseElasticV2::~G4DiffuseElasticV2 ( )
virtual

Definition at line 115 of file G4DiffuseElasticV2.cc.

116{
117 if ( fEnergyVector )
118 {
119 delete fEnergyVector;
120 fEnergyVector = 0;
121 }
122}

Member Function Documentation

◆ BesselJone()

G4double G4DiffuseElasticV2::BesselJone ( G4double z)
inline

Definition at line 267 of file G4DiffuseElasticV2.hh.

268{
269 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
270
271 modvalue = std::fabs(value);
272
273 if ( modvalue < 8.0 )
274 {
275 value2 = value*value;
276
277 fact1 = value*(72362614232.0 + value2*(-7895059235.0
278 + value2*( 242396853.1
279 + value2*(-2972611.439
280 + value2*( 15704.48260
281 + value2*(-30.16036606 ) ) ) ) ) );
282
283 fact2 = 144725228442.0 + value2*(2300535178.0
284 + value2*(18583304.74
285 + value2*(99447.43394
286 + value2*(376.9991397
287 + value2*1.0 ) ) ) );
288 bessel = fact1/fact2;
289 }
290 else
291 {
292 arg = 8.0/modvalue;
293
294 value2 = arg*arg;
295
296 shift = modvalue - 2.356194491;
297
298 fact1 = 1.0 + value2*( 0.183105e-2
299 + value2*(-0.3516396496e-4
300 + value2*(0.2457520174e-5
301 + value2*(-0.240337019e-6 ) ) ) );
302
303 fact2 = 0.04687499995 + value2*(-0.2002690873e-3
304 + value2*( 0.8449199096e-5
305 + value2*(-0.88228987e-6
306 + value2*0.105787412e-6 ) ) );
307
308 bessel = std::sqrt( 0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2);
309
310 if (value < 0.0) bessel = -bessel;
311 }
312 return bessel;
313}
double G4double
Definition G4Types.hh:83

Referenced by BesselOneByArg(), and GetDiffElasticSumProbA().

◆ BesselJzero()

G4double G4DiffuseElasticV2::BesselJzero ( G4double z)
inline

Definition at line 215 of file G4DiffuseElasticV2.hh.

216{
217 G4double modvalue, value2, fact1, fact2, arg, shift, bessel;
218
219 modvalue = std::fabs(value);
220
221 if ( value < 8.0 && value > -8.0 )
222 {
223 value2 = value*value;
224
225 fact1 = 57568490574.0 + value2*(-13362590354.0
226 + value2*( 651619640.7
227 + value2*(-11214424.18
228 + value2*( 77392.33017
229 + value2*(-184.9052456 ) ) ) ) );
230
231 fact2 = 57568490411.0 + value2*( 1029532985.0
232 + value2*( 9494680.718
233 + value2*(59272.64853
234 + value2*(267.8532712
235 + value2*1.0 ) ) ) );
236
237 bessel = fact1/fact2;
238 }
239 else
240 {
241 arg = 8.0/modvalue;
242
243 value2 = arg*arg;
244
245 shift = modvalue-0.785398164;
246
247 fact1 = 1.0 + value2*(-0.1098628627e-2
248 + value2*(0.2734510407e-4
249 + value2*(-0.2073370639e-5
250 + value2*0.2093887211e-6 ) ) );
251
252 fact2 = -0.1562499995e-1 + value2*(0.1430488765e-3
253 + value2*(-0.6911147651e-5
254 + value2*(0.7621095161e-6
255 - value2*0.934945152e-7 ) ) );
256
257 bessel = std::sqrt(0.636619772/modvalue)*(std::cos(shift)*fact1 - arg*std::sin(shift)*fact2 );
258 }
259 return bessel;
260}

Referenced by GetDiffElasticSumProbA().

◆ BesselOneByArg()

G4double G4DiffuseElasticV2::BesselOneByArg ( G4double z)
inline

Definition at line 342 of file G4DiffuseElasticV2.hh.

343{
344 G4double x2, result;
345
346 if( std::fabs(x) < 0.01 )
347 {
348 x *= 0.5;
349 x2 = x*x;
350 result = 2. - x2 + x2*x2/6.;
351 }
352 else
353 {
354 result = BesselJone(x)/x;
355 }
356 return result;
357}
G4double BesselJone(G4double z)

Referenced by GetDiffElasticSumProbA().

◆ BuildAngleTable()

void G4DiffuseElasticV2::BuildAngleTable ( )

Definition at line 435 of file G4DiffuseElasticV2.cc.

436{
437 G4double partMom, kinE, a = 0., z = fParticle->GetPDGCharge(), m1 = fParticle->GetPDGMass();
438 G4double alpha1, alpha2, alphaMax, alphaCoulomb, delta = 0., sum = 0.;
439
441
442 fEnergyAngleVector = new std::vector<std::vector<G4double>*>;
443 fEnergySumVector = new std::vector<std::vector<G4double>*>;
444
445 for( G4int i = 0; i < fEnergyBin; ++i)
446 {
447 kinE = fEnergyVector->Energy(i);
448 partMom = std::sqrt( kinE*(kinE + 2*m1) );
449
450 fWaveVector = partMom/hbarc;
451
452 G4double kR = fWaveVector*fNuclearRadius;
453 G4double kRmax = 18.6; // 10.6; 10.6, 18, 10.174; ~ 3 maxima of J1 or 15., 25.
454 G4double kRcoul = 1.9; // 1.2; 1.4, 2.5; // on the first slope of J1
455
456 alphaMax = kRmax/kR;
457
458 if ( alphaMax >= CLHEP::pi ) alphaMax = CLHEP::pi; // vmg21.10.15
459
460 alphaCoulomb = kRcoul/kR;
461
462 if( z )
463 {
464 a = partMom/m1; // beta*gamma for m1
465 fBeta = a/std::sqrt(1+a*a);
466 fZommerfeld = CalculateZommerfeld( fBeta, z, fAtomicNumber);
467 fAm = CalculateAm( partMom, fZommerfeld, fAtomicNumber);
468 fAddCoulomb = true;
469 }
470
471 std::vector<G4double>* angleVector = new std::vector<G4double>(fAngleBin);
472 std::vector<G4double>* sumVector = new std::vector<G4double>(fAngleBin);
473
474
475 G4double delth = alphaMax/fAngleBin;
476
477 sum = 0.;
478
479 for(G4int j = (G4int)fAngleBin-1; j >= 0; --j)
480 {
481 alpha1 = delth*j;
482 alpha2 = alpha1 + delth;
483
484 if( fAddCoulomb && ( alpha2 < alphaCoulomb)) fAddCoulomb = false;
485
486 delta = integral.Legendre10(this, &G4DiffuseElasticV2::GetIntegrandFunction, alpha1, alpha2);
487
488 sum += delta;
489
490 (*angleVector)[j] = alpha1;
491 (*sumVector)[j] = sum;
492
493 }
494 fEnergyAngleVector->push_back(angleVector);
495 fEnergySumVector->push_back(sumVector);
496
497 }
498 return;
499}
int G4int
Definition G4Types.hh:85
const G4double alpha2
G4double CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2)
G4double CalculateAm(G4double momentum, G4double n, G4double Z)
G4double GetIntegrandFunction(G4double theta)
G4double Energy(const std::size_t index) const

Referenced by Initialise(), and InitialiseOnFly().

◆ CalculateAm()

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

Definition at line 375 of file G4DiffuseElasticV2.hh.

376{
377 G4double k = momentum/CLHEP::hbarc;
378 G4double ch = 1.13 + 3.76*n*n;
379 G4double zn = 1.77*k*(1.0/G4Pow::GetInstance()->A13(Z))*CLHEP::Bohr_radius;
380 G4double zn2 = zn*zn;
381 fAm = ch/zn2;
382
383 return fAm;
384}
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double A13(G4double A) const
Definition G4Pow.cc:116

Referenced by BuildAngleTable().

◆ CalculateNuclearRad()

G4double G4DiffuseElasticV2::CalculateNuclearRad ( G4double A)
inline

Definition at line 390 of file G4DiffuseElasticV2.hh.

391{
392 G4double R, r0, a11, a12, a13, a2, a3;
393
394 a11 = 1.26; // 1.08, 1.16
395 a12 = 1.; // 1.08, 1.16
396 a13 = 1.12; // 1.08, 1.16
397 a2 = 1.1;
398 a3 = 1.;
399
400 // Special rms radii for light nucleii
401
402 if (A < 50.)
403 {
404 if (std::abs(A-1.) < 0.5) return 0.89*CLHEP::fermi; // p
405 else if(std::abs(A-2.) < 0.5) return 2.13*CLHEP::fermi; // d
406 else if( // std::abs(Z-1.) < 0.5 &&
407 std::abs(A-3.) < 0.5) return 1.80*CLHEP::fermi; // t
408
409 // else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96CLHEP::fermi; // He3
410 else if( // std::abs(Z-2.) < 0.5 &&
411 std::abs(A-4.) < 0.5) return 1.68*CLHEP::fermi; // He4
412
413 else if( // std::abs(Z-3.) < 0.5
414 std::abs(A-7.) < 0.5 ) return 2.40*CLHEP::fermi; // Li7
415 else if( // std::abs(Z-4.) < 0.5
416 std::abs(A-9.) < 0.5) return 2.51*CLHEP::fermi; // Be9
417
418 else if( 10. < A && A <= 16. ) r0 = a11*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi; // 1.08CLHEP::fermi;
419 else if( 15. < A && A <= 20. ) r0 = a12*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
420 else if( 20. < A && A <= 30. ) r0 = a13*( 1 - (1.0/G4Pow::GetInstance()->A23(A)) )*CLHEP::fermi;
421 else r0 = a2*CLHEP::fermi;
422
423 R = r0*G4Pow::GetInstance()->A13(A);
424 }
425 else
426 {
427 r0 = a3*CLHEP::fermi;
428
429 R = r0*G4Pow::GetInstance()->powA(A, 0.27);
430 }
431 fNuclearRadius = R;
432
433 return R;
434}
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 Initialise(), and InitialiseOnFly().

◆ CalculateZommerfeld()

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

Definition at line 364 of file G4DiffuseElasticV2.hh.

365{
366 fZommerfeld = CLHEP::fine_structure_const*Z1*Z2/beta;
367
368 return fZommerfeld;
369}

Referenced by BuildAngleTable().

◆ DampFactor()

G4double G4DiffuseElasticV2::DampFactor ( G4double z)
inline

Definition at line 319 of file G4DiffuseElasticV2.hh.

320{
321 G4double df;
322 G4double f2 = 2., f3 = 6., f4 = 24.; // first factorials
323
324 // x *= pi;
325
326 if( std::fabs(x) < 0.01 )
327 {
328 df = 1./(1. + x/f2 + x*x/f3 + x*x*x/f4);
329 }
330 else
331 {
332 df = x/std::sinh(x);
333 }
334 return df;
335}

Referenced by GetDiffElasticSumProbA().

◆ GetDiffElasticSumProbA()

G4double G4DiffuseElasticV2::GetDiffElasticSumProbA ( G4double alpha)

Definition at line 164 of file G4DiffuseElasticV2.cc.

165{
166
167 G4double sigma, bzero, bzero2, bonebyarg, bonebyarg2, damp, damp2;
168 G4double delta, diffuse, gamma;
169 G4double e1, e2, bone, bone2;
170
171 // G4double wavek = momentum/hbarc; // wave vector
172 // G4double r0 = 1.08*fermi;
173 // G4double rad = r0*G4Pow::GetInstance()->A13(A);
174
175 G4double kr = fWaveVector*fNuclearRadius; // wavek*rad;
176 G4double kr2 = kr*kr;
177 G4double krt = kr*theta;
178
179 bzero = BesselJzero(krt);
180 bzero2 = bzero*bzero;
181 bone = BesselJone(krt);
182 bone2 = bone*bone;
183 bonebyarg = BesselOneByArg(krt);
184 bonebyarg2 = bonebyarg*bonebyarg;
185
186 if ( fParticle == theProton )
187 {
188 diffuse = 0.63*fermi;
189 gamma = 0.3*fermi;
190 delta = 0.1*fermi*fermi;
191 e1 = 0.3*fermi;
192 e2 = 0.35*fermi;
193 }
194 else if ( fParticle == theNeutron )
195 {
196 diffuse = 0.63*fermi;
197 gamma = 0.3*fermi;
198 delta = 0.1*fermi*fermi;
199 e1 = 0.3*fermi;
200 e2 = 0.35*fermi;
201 }
202 else // as proton, if were not defined
203 {
204 diffuse = 0.63*fermi;
205 gamma = 0.3*fermi;
206 delta = 0.1*fermi*fermi;
207 e1 = 0.3*fermi;
208 e2 = 0.35*fermi;
209 }
210
211 G4double lambda = 15; // 15 ok
212 // G4double kgamma = fWaveVector*gamma; // wavek*delta;
213 G4double kgamma = lambda*(1.-G4Exp(-fWaveVector*gamma/lambda)); // wavek*delta;
214
215 if( fAddCoulomb ) // add Coulomb correction
216 {
217 G4double sinHalfTheta = std::sin(0.5*theta);
218 G4double sinHalfTheta2 = sinHalfTheta*sinHalfTheta;
219
220 kgamma += 0.5*fZommerfeld/kr/(sinHalfTheta2+fAm); // correction at J0()
221 }
222
223 G4double kgamma2 = kgamma*kgamma;
224
225 // G4double dk2t = delta*fWaveVector*fWaveVector*theta; // delta*wavek*wavek*theta;
226 // G4double dk2t2 = dk2t*dk2t;
227
228 // G4double pikdt = pi*fWaveVector*diffuse*theta;// pi*wavek*diffuse*theta;
229 G4double pikdt = lambda*(1. - G4Exp( -pi*fWaveVector*diffuse*theta/lambda ) ); // wavek*delta;
230
231 damp = DampFactor( pikdt );
232 damp2 = damp*damp;
233
234 G4double mode2k2 = ( e1*e1 + e2*e2 )*fWaveVector*fWaveVector;
235 G4double e2dk3t = -2.*e2*delta*fWaveVector*fWaveVector*fWaveVector*theta;
236
237 sigma = kgamma2;
238 // sigma += dk2t2;
239 sigma *= bzero2;
240 sigma += mode2k2*bone2;
241 sigma += e2dk3t*bzero*bone;
242
243 // sigma += kr2*(1 + 8.*fZommerfeld*fZommerfeld/kr2)*bonebyarg2; // correction at J1()/()
244 sigma += kr2*bonebyarg2; // correction at J1()/()
245
246 sigma *= damp2; // *rad*rad;
247
248 return sigma;
249}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double BesselJzero(G4double z)
G4double DampFactor(G4double z)
G4double BesselOneByArg(G4double z)

Referenced by GetIntegrandFunction().

◆ GetIntegrandFunction()

G4double G4DiffuseElasticV2::GetIntegrandFunction ( G4double theta)

Definition at line 257 of file G4DiffuseElasticV2.cc.

258{
259 G4double result;
260
261 result = GetDiffElasticSumProbA(alpha) * 2 * CLHEP::pi * std::sin(alpha);
262
263 return result;
264}
G4double GetDiffElasticSumProbA(G4double alpha)

Referenced by BuildAngleTable().

◆ GetNuclearRadius()

G4double G4DiffuseElasticV2::GetNuclearRadius ( )
inline

Definition at line 129 of file G4DiffuseElasticV2.hh.

129{return fNuclearRadius;};

◆ GetScatteringAngle()

G4double G4DiffuseElasticV2::GetScatteringAngle ( G4int iMomentum,
unsigned long iAngle,
G4double position )

Definition at line 506 of file G4DiffuseElasticV2.cc.

507{
508 G4double x1, x2, y1, y2, randAngle = 0;
509
510 if( iAngle == 0 )
511 {
512 randAngle = (*(*fEnergyAngleVector)[iMomentum])[iAngle];
513 }
514 else
515 {
516 if ( iAngle >= (*fEnergyAngleVector)[iMomentum]->size() )
517 {
518 iAngle = (*fEnergyAngleVector)[iMomentum]->size() - 1;
519 }
520
521 y1 = (*(*fEnergySumVector)[iMomentum])[iAngle-1];
522 y2 = (*(*fEnergySumVector)[iMomentum])[iAngle];
523
524 x1 = (*(*fEnergyAngleVector)[iMomentum])[iAngle-1];
525 x2 = (*(*fEnergyAngleVector)[iMomentum])[iAngle];
526
527 if ( x1 == x2 ) randAngle = x2;
528 else
529 {
530 if ( y1 == y2 ) randAngle = x1 + ( x2 - x1 )*G4UniformRand();
531 else
532 {
533 randAngle = x1 + ( position - y1 )*( x2 - x1 )/( y2 - y1 );
534 }
535 }
536 }
537
538 return randAngle;
539}
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by SampleTableThetaCMS().

◆ Initialise()

void G4DiffuseElasticV2::Initialise ( )

Definition at line 128 of file G4DiffuseElasticV2.cc.

129{
130
131 const G4ElementTable* theElementTable = G4Element::GetElementTable();
132
133 std::size_t jEl, numOfEl = G4Element::GetNumberOfElements();
134
135 for( jEl = 0; jEl < numOfEl; ++jEl) // application element loop
136 {
137 fAtomicNumber = (*theElementTable)[jEl]->GetZ(); // atomic number
138 fAtomicWeight = G4NistManager::Instance()->GetAtomicMassAmu( static_cast< G4int >( fAtomicNumber ) );
139 fNuclearRadius = CalculateNuclearRad(fAtomicWeight);
140
141 if( verboseLevel > 0 )
142 {
143 G4cout<<"G4DiffuseElasticV2::Initialise() the element: "
144 <<(*theElementTable)[jEl]->GetName()<<G4endl;
145 }
146 fElementNumberVector.push_back(fAtomicNumber);
147 fElementNameVector.push_back((*theElementTable)[jEl]->GetName());
148
150
151 fEnergyAngleVectorBank.push_back(fEnergyAngleVector);
152 fEnergySumVectorBank.push_back(fEnergySumVector);
153
154 }
155 return;
156}
std::vector< G4Element * > G4ElementTable
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double CalculateNuclearRad(G4double A)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static size_t GetNumberOfElements()
Definition G4Element.cc:393
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const

◆ InitialiseOnFly()

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

Definition at line 408 of file G4DiffuseElasticV2.cc.

409{
410 fAtomicNumber = Z; // atomic number
411 fAtomicWeight = G4NistManager::Instance()->GetAtomicMassAmu( static_cast< G4int >( Z ) );
412
413 fNuclearRadius = CalculateNuclearRad(fAtomicWeight);
414
415 if( verboseLevel > 0 )
416 {
417 G4cout<<"G4DiffuseElasticV2::InitialiseOnFly() the element with Z = "
418 <<Z<<"; and A = "<<A<<G4endl;
419 }
420 fElementNumberVector.push_back(fAtomicNumber);
421
423
424 fEnergyAngleVectorBank.push_back(fEnergyAngleVector);
425 fEnergySumVectorBank.push_back(fEnergySumVector);
426
427 return;
428}

Referenced by SampleTableThetaCMS().

◆ IsApplicable()

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

Reimplemented from G4HadronicInteraction.

Definition at line 170 of file G4DiffuseElasticV2.hh.

172{
173 if( ( projectile.GetDefinition() == G4Proton::Proton() ||
174 projectile.GetDefinition() == G4Neutron::Neutron() ||
175 projectile.GetDefinition() == G4PionPlus::PionPlus() ||
176 projectile.GetDefinition() == G4PionMinus::PionMinus() ||
177 projectile.GetDefinition() == G4KaonPlus::KaonPlus() ||
178 projectile.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
179
180 nucleus.GetZ_asInt() >= 2 ) return true;
181 else return false;
182}
const G4ParticleDefinition * GetDefinition() const
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93

◆ NeutronTuniform()

G4double G4DiffuseElasticV2::NeutronTuniform ( G4int Z)

Definition at line 310 of file G4DiffuseElasticV2.cc.

311{
312 G4double elZ = G4double(Z);
313 elZ -= 1.;
314 G4double Tkin = 12.*G4Exp(-elZ/10.) + 1.;
315
316 return Tkin;
317}

Referenced by SampleInvariantT().

◆ SampleInvariantT()

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

Reimplemented from G4HadronicInteraction.

Definition at line 273 of file G4DiffuseElasticV2.cc.

275{
276 fParticle = aParticle;
277 G4double m1 = fParticle->GetPDGMass(), t;
278 G4double totElab = std::sqrt(m1*m1+p*p);
280 G4LorentzVector lv1(p,0.0,0.0,totElab);
281 G4LorentzVector lv(0.0,0.0,0.0,mass2);
282 lv += lv1;
283
284 G4ThreeVector bst = lv.boostVector();
285 lv1.boost(-bst);
286
287 G4ThreeVector p1 = lv1.vect();
288 G4double momentumCMS = p1.mag();
289
290 if( aParticle == theNeutron)
291 {
292 G4double Tmax = NeutronTuniform( Z );
293 G4double pCMS2 = momentumCMS*momentumCMS;
294 G4double Tkin = std::sqrt(pCMS2+m1*m1)-m1;
295
296 if( Tkin <= Tmax )
297 {
298 t = 4.*pCMS2*G4UniformRand();
299 return t;
300 }
301 }
302
303 t = SampleTableT( aParticle, momentumCMS, G4double(Z), G4double(A) ); // sample theta in cms
304
305 return t;
306}
double mag() const
G4double SampleTableT(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)
G4double NeutronTuniform(G4int Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ SampleTableT()

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

Definition at line 324 of file G4DiffuseElasticV2.cc.

326{
327 G4double alpha = SampleTableThetaCMS( aParticle, p, Z, A); // sample theta in cms
328 G4double t = 2*p*p*( 1 - std::cos(alpha) ); // -t !!!
329
330 return t;
331}
G4double SampleTableThetaCMS(const G4ParticleDefinition *aParticle, G4double p, G4double Z, G4double A)

Referenced by SampleInvariantT().

◆ SampleTableThetaCMS()

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

Definition at line 339 of file G4DiffuseElasticV2.cc.

341{
342 std::size_t iElement;
343 G4int iMomentum;
344 unsigned long iAngle = 0;
345 G4double randAngle, position, theta1, theta2, E1, E2, W1, W2, W;
346 G4double m1 = particle->GetPDGMass();
347
348 for(iElement = 0; iElement < fElementNumberVector.size(); iElement++)
349 {
350 if( std::fabs(Z - fElementNumberVector[iElement]) < 0.5) break;
351 }
352
353 if ( iElement == fElementNumberVector.size() )
354 {
355 InitialiseOnFly(Z,A); // table preparation, if needed
356 }
357
358 fEnergyAngleVector = fEnergyAngleVectorBank[iElement];
359 fEnergySumVector = fEnergySumVectorBank[iElement];
360
361
362 G4double kinE = std::sqrt(momentum*momentum + m1*m1) - m1;
363
364 iMomentum = G4int(fEnergyVector->FindBin(kinE,1000) + 1);
365
366 position = (*(*fEnergySumVector)[iMomentum])[0]*G4UniformRand();
367
368 for(iAngle = 0; iAngle < fAngleBin; ++iAngle)
369 {
370 if (position > (*(*fEnergySumVector)[iMomentum])[iAngle]) break;
371 }
372
373
374 if (iMomentum == fEnergyBin -1 || iMomentum == 0 ) // the table edges
375 {
376 randAngle = GetScatteringAngle(iMomentum, iAngle, position);
377 }
378 else // kinE inside between energy table edges
379 {
380 theta2 = GetScatteringAngle(iMomentum, iAngle, position);
381
382 E2 = fEnergyVector->Energy(iMomentum);
383
384 iMomentum--;
385
386 theta1 = GetScatteringAngle(iMomentum, iAngle, position);
387
388 E1 = fEnergyVector->Energy(iMomentum);
389
390 W = 1.0/(E2 - E1);
391 W1 = (E2 - kinE)*W;
392 W2 = (kinE - E1)*W;
393
394 randAngle = W1*theta1 + W2*theta2;
395 }
396
397
398
399 if(randAngle < 0.) randAngle = 0.;
400
401 return randAngle;
402}
G4double GetScatteringAngle(G4int iMomentum, unsigned long iAngle, G4double position)
void InitialiseOnFly(G4double Z, G4double A)
std::size_t FindBin(const G4double energy, std::size_t idx) const
#define W
Definition crc32.c:85

Referenced by SampleTableT().

◆ SampleThetaCMS()

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

◆ SampleThetaLab()

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

◆ SetHEModelLowLimit()

void G4DiffuseElasticV2::SetHEModelLowLimit ( G4double value)
inline

Definition at line 194 of file G4DiffuseElasticV2.hh.

195{
196 lowEnergyLimitHE = value;
197}

◆ SetLowestEnergyLimit()

void G4DiffuseElasticV2::SetLowestEnergyLimit ( G4double value)
inline

Definition at line 204 of file G4DiffuseElasticV2.hh.

205{
206 lowestEnergyLimit = value;
207}

◆ SetPlabLowLimit()

void G4DiffuseElasticV2::SetPlabLowLimit ( G4double value)
inline

Definition at line 189 of file G4DiffuseElasticV2.hh.

190{
191 plabLowLimit = value;
192}

◆ SetQModelLowLimit()

void G4DiffuseElasticV2::SetQModelLowLimit ( G4double value)
inline

Definition at line 199 of file G4DiffuseElasticV2.hh.

200{
201 lowEnergyLimitQ = value;
202}

◆ SetRecoilKinEnergyLimit()

void G4DiffuseElasticV2::SetRecoilKinEnergyLimit ( G4double value)
inline

Definition at line 184 of file G4DiffuseElasticV2.hh.

185{
186 lowEnergyRecoilLimit = value;
187}

◆ ThetaCMStoThetaLab()

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

Definition at line 551 of file G4DiffuseElasticV2.cc.

553{
554 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
555 G4double m1 = theParticle->GetPDGMass();
556 G4LorentzVector lv1 = aParticle->Get4Momentum();
557 G4LorentzVector lv(0.0,0.0,0.0,tmass);
558
559 lv += lv1;
560
561 G4ThreeVector bst = lv.boostVector();
562
563 lv1.boost(-bst);
564
565 G4ThreeVector p1 = lv1.vect();
566 G4double ptot = p1.mag();
567
568 G4double phi = G4UniformRand()*twopi;
569 G4double cost = std::cos(thetaCMS);
570 G4double sint;
571
572 if( cost >= 1.0 )
573 {
574 cost = 1.0;
575 sint = 0.0;
576 }
577 else if( cost <= -1.0)
578 {
579 cost = -1.0;
580 sint = 0.0;
581 }
582 else
583 {
584 sint = std::sqrt((1.0-cost)*(1.0+cost));
585 }
586 if (verboseLevel>1)
587 {
588 G4cout << "cos(tcms)=" << cost << " std::sin(tcms)=" << sint << G4endl;
589 }
590 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
591 v1 *= ptot;
592 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(ptot*ptot + m1*m1));
593
594 nlv1.boost(bst);
595
596 G4ThreeVector np1 = nlv1.vect();
597
598 G4double thetaLab = np1.theta();
599
600 return thetaLab;
601}
double theta() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const

◆ ThetaLabToThetaCMS()

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

Definition at line 609 of file G4DiffuseElasticV2.cc.

611{
612 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
613 G4double m1 = theParticle->GetPDGMass();
614 G4double plab = aParticle->GetTotalMomentum();
615 G4LorentzVector lv1 = aParticle->Get4Momentum();
616 G4LorentzVector lv(0.0,0.0,0.0,tmass);
617
618 lv += lv1;
619
620 G4ThreeVector bst = lv.boostVector();
621
622 G4double phi = G4UniformRand()*twopi;
623 G4double cost = std::cos(thetaLab);
624 G4double sint;
625
626 if( cost >= 1.0 )
627 {
628 cost = 1.0;
629 sint = 0.0;
630 }
631 else if( cost <= -1.0)
632 {
633 cost = -1.0;
634 sint = 0.0;
635 }
636 else
637 {
638 sint = std::sqrt((1.0-cost)*(1.0+cost));
639 }
640 if (verboseLevel>1)
641 {
642 G4cout << "cos(tlab)=" << cost << " std::sin(tlab)=" << sint << G4endl;
643 }
644 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
645 v1 *= plab;
646 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),std::sqrt(plab*plab + m1*m1));
647
648 nlv1.boost(-bst);
649
650 G4ThreeVector np1 = nlv1.vect();
651 G4double thetaCMS = np1.theta();
652
653 return thetaCMS;
654}
G4double GetTotalMomentum() const

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