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

#include <G4ParticleHPThermalScattering.hh>

Inheritance diagram for G4ParticleHPThermalScattering:

Public Member Functions

 G4ParticleHPThermalScattering ()
 ~G4ParticleHPThermalScattering () override
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus) override
const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const override
void AddUserThermalScatteringFile (const G4String &, const G4String &)
void BuildPhysicsTable (const G4ParticleDefinition &) override
void ModelDescription (std::ostream &outFile) const override
Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
virtual ~G4HadronicInteraction ()
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 std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
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 G4HadronicInteraction
G4HadFinalState theParticleChange
G4int verboseLevel
G4double theMinEnergy
G4double theMaxEnergy
G4bool isBlocked

Detailed Description

Definition at line 83 of file G4ParticleHPThermalScattering.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPThermalScattering()

G4ParticleHPThermalScattering::G4ParticleHPThermalScattering ( )

Definition at line 56 of file G4ParticleHPThermalScattering.cc.

57 : G4HadronicInteraction("NeutronHPThermalScattering")
58{
59 theHPElastic = new G4ParticleHPElastic();
60
61 SetMinEnergy(0. * eV);
62 SetMaxEnergy(4 * eV);
63 theXSection = new G4ParticleHPThermalScatteringData();
64
65 nMaterial = 0;
66 nElement = 0;
67}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4ParticleHPThermalScattering()

G4ParticleHPThermalScattering::~G4ParticleHPThermalScattering ( )
override

Definition at line 69 of file G4ParticleHPThermalScattering.cc.

70{
71 delete theHPElastic;
72}

Member Function Documentation

◆ AddUserThermalScatteringFile()

void G4ParticleHPThermalScattering::AddUserThermalScatteringFile ( const G4String & nameG4Element,
const G4String & filename )

Definition at line 1188 of file G4ParticleHPThermalScattering.cc.

1190{
1191 names.AddThermalElement(nameG4Element, filename);
1192 theXSection->AddUserThermalScatteringFile(nameG4Element, filename);
1193 buildPhysicsTable();
1194}

◆ ApplyYourself()

G4HadFinalState * G4ParticleHPThermalScattering::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & aTargetNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 305 of file G4ParticleHPThermalScattering.cc.

307{
308 // Select Element > Reaction >
309
310 const G4Material* theMaterial = aTrack.GetMaterial();
311 G4double aTemp = theMaterial->GetTemperature();
312 auto n = (G4int)theMaterial->GetNumberOfElements();
313
314 G4bool findThermalElement = false;
315 G4int ielement;
316 const G4Element* theElement = nullptr;
317 for (G4int i = 0; i < n; ++i) {
318 theElement = theMaterial->GetElement(i);
319 // Select target element
320 if (aNucleus.GetZ_asInt() == (G4int)(theElement->GetZ() + 0.5)) {
321 // Check Applicability of Thermal Scattering
322 if (getTS_ID(nullptr, theElement) != -1) {
323 ielement = getTS_ID(nullptr, theElement);
324 findThermalElement = true;
325 break;
326 }
327 if (getTS_ID(theMaterial, theElement) != -1) {
328 ielement = getTS_ID(theMaterial, theElement);
329 findThermalElement = true;
330 break;
331 }
332 }
333 }
334
335 if (findThermalElement) {
336 // Select Reaction (Inelastic, coherent, incoherent)
337 const G4ParticleDefinition* pd = aTrack.GetDefinition();
338 auto dp = new G4DynamicParticle(pd, aTrack.Get4Momentum());
339 G4double total = theXSection->GetCrossSection(dp, theElement, theMaterial);
340 G4double inelastic = theXSection->GetInelasticCrossSection(dp, theElement, theMaterial);
341
342 auto inelELM = inelasticFSs->find(ielement);
343 G4bool mayBeInel = (inelELM != inelasticFSs->end());
344 auto coheELM = coherentFSs->find(ielement);
345 G4bool mayBeCohe = (coheELM != coherentFSs->end());
346 auto incoELM = incoherentFSs->find(ielement);
347 G4bool mayBeInco = (incoELM != incoherentFSs->end());
348
349 G4double random = G4UniformRand();
350 if ((total*random <= inelastic && mayBeInel) || (!mayBeCohe && !mayBeInco)) {
351 // Inelastic
352
353 std::vector<G4double> v_temp;
354 v_temp.clear();
355 for (auto it = inelELM->second->begin(); it != inelELM->second->end(); ++it)
356 {
357 v_temp.push_back(it->first);
358 }
359
360 std::pair<G4double, G4double> tempLH = find_LH(aTemp, &v_temp);
361 //
362 // For T_L aNEP_EPM_TL and T_H aNEP_EPM_TH
363 //
364 std::vector<E_P_E_isoAng*>* vNEP_EPM_TL = nullptr;
365 std::vector<E_P_E_isoAng*>* vNEP_EPM_TH = nullptr;
366
367 if (tempLH.first != 0.0 && tempLH.second != 0.0) {
368 vNEP_EPM_TL = inelELM->second->find(tempLH.first / kelvin)->second;
369 vNEP_EPM_TH = inelELM->second->find(tempLH.second / kelvin)->second;
370 }
371 else if (tempLH.first == 0.0) {
372 auto itm = inelELM->second->begin();
373 vNEP_EPM_TL = itm->second;
374 ++itm;
375 vNEP_EPM_TH = itm->second;
376 tempLH.first = tempLH.second;
377 tempLH.second = itm->first;
378 }
379 else if (tempLH.second == 0.0) {
380 auto itm = inelELM->second->end();
381 --itm;
382 vNEP_EPM_TH = itm->second;
383 --itm;
384 vNEP_EPM_TL = itm->second;
385 tempLH.second = tempLH.first;
386 tempLH.first = itm->first;
387 }
388
389 G4double sE = 0., mu = 1.0;
390
391 // New Geant4 method - Stochastic temperature interpolation of the final state
392 // (continuous temperature interpolation was used previously)
393 std::pair<G4double, G4double> secondaryParam;
394 G4double rand_temp = G4UniformRand();
395 if (rand_temp < (aTemp - tempLH.first) / (tempLH.second - tempLH.first))
396 secondaryParam = sample_inelastic_E_mu(aTrack.GetKineticEnergy(), vNEP_EPM_TH);
397 else
398 secondaryParam = sample_inelastic_E_mu(aTrack.GetKineticEnergy(), vNEP_EPM_TL);
399
400 sE = secondaryParam.first;
401 mu = secondaryParam.second;
402
403 // set
404 theParticleChange.SetEnergyChange(sE);
405 G4double phi = CLHEP::twopi * G4UniformRand();
406 G4double sint = std::sqrt(1 - mu * mu);
407 theParticleChange.SetMomentumChange(sint * std::cos(phi), sint * std::sin(phi), mu);
408 }
409 else if (mayBeCohe &&
410 (!mayBeInco || random*total
411 <= (inelastic + theXSection->GetCoherentCrossSection(dp, theElement, theMaterial))))
412 {
413 // Coherent Elastic
414
415 G4double E = aTrack.GetKineticEnergy();
416
417 // T_L and T_H
418 std::vector<G4double> v_temp;
419 v_temp.clear();
420 for (auto it = coheELM->second->begin(); it != coheELM->second->end(); ++it)
421 {
422 v_temp.push_back(it->first);
423 }
424
425 // T_L T_H
426 std::pair<G4double, G4double> tempLH = find_LH(aTemp, &v_temp);
427 //
428 //
429 // For T_L anEPM_TL and T_H anEPM_TH
430 //
431 std::vector<std::pair<G4double, G4double>*>* pvE_p_TL = nullptr;
432 std::vector<std::pair<G4double, G4double>*>* pvE_p_TH = nullptr;
433
434 if (tempLH.first != 0.0 && tempLH.second != 0.0) {
435 pvE_p_TL = coheELM->second->find(tempLH.first / kelvin)->second;
436 pvE_p_TH = coheELM->second->find(tempLH.first / kelvin)->second;
437 }
438 else if (tempLH.first == 0.0) {
439 pvE_p_TL = coheELM->second->find(v_temp[0])->second;
440 pvE_p_TH = coheELM->second->find(v_temp[1])->second;
441 tempLH.first = tempLH.second;
442 tempLH.second = v_temp[1];
443 }
444 else if (tempLH.second == 0.0) {
445 pvE_p_TH = coheELM->second->find(v_temp.back())->second;
446 auto itv = v_temp.cend();
447 --itv;
448 --itv;
449 pvE_p_TL = coheELM->second->find(*itv)->second;
450 tempLH.second = tempLH.first;
451 tempLH.first = *itv;
452 }
453 else {
454 // tempLH.first == 0.0 && tempLH.second
455 throw G4HadronicException(
456 __FILE__, __LINE__,
457 "A problem is found in Thermal Scattering Data! Unexpected temperature values in data");
458 }
459
460 std::vector<G4double> vE_T;
461 std::vector<G4double> vp_T;
462
463 auto n1 = (G4int)pvE_p_TL->size();
464
465 // New Geant4 method - Stochastic interpolation of the final state
466 std::vector<std::pair<G4double, G4double>*>* pvE_p_T_sampled;
467 G4double rand_temp = G4UniformRand();
468 if (rand_temp < (aTemp - tempLH.first) / (tempLH.second - tempLH.first))
469 pvE_p_T_sampled = pvE_p_TH;
470 else
471 pvE_p_T_sampled = pvE_p_TL;
472
473 // 171005 fix bug, contribution from H.N. TRAN@CEA
474 for (G4int i = 0; i < n1; ++i) {
475 vE_T.push_back((*pvE_p_T_sampled)[i]->first);
476 vp_T.push_back((*pvE_p_T_sampled)[i]->second);
477 }
478
479 G4int j = 0;
480 for (G4int i = 1; i < n1; ++i) {
481 if (E / eV < vE_T[i]) {
482 j = i - 1;
483 break;
484 }
485 }
486
487 G4double rand_for_mu = G4UniformRand();
488
489 G4int k = 0;
490 for (G4int i = 0; i <= j; ++i) {
491 G4double Pi = vp_T[i] / vp_T[j];
492 if (rand_for_mu < Pi) {
493 k = i;
494 break;
495 }
496 }
497
498 G4double Ei = vE_T[k];
499
500 G4double mu = 1 - 2 * Ei / (E / eV);
501
502 if (mu < -1.0) mu = -1.0;
503
504 theParticleChange.SetEnergyChange(E);
505 G4double phi = CLHEP::twopi * G4UniformRand();
506 G4double sint = std::sqrt(1 - mu * mu);
507 theParticleChange.SetMomentumChange(sint * std::cos(phi), sint * std::sin(phi), mu);
508 }
509 else if (mayBeInco) {
510 // InCoherent Elastic
511
512 // T_L and T_H
513 std::vector<G4double> v_temp;
514 v_temp.clear();
515 for (auto it = incoELM->second->cbegin(); it != incoELM->second->cend(); ++it)
516 {
517 v_temp.push_back(it->first);
518 }
519
520 // T_L T_H
521 std::pair<G4double, G4double> tempLH = find_LH(aTemp, &v_temp);
522
523 //
524 // For T_L anEPM_TL and T_H anEPM_TH
525 //
526
527 E_isoAng anEPM_TL_E;
528 E_isoAng anEPM_TH_E;
529
530 if (tempLH.first != 0.0 && tempLH.second != 0.0) {
531 // Interpolate TL and TH
532 anEPM_TL_E = create_E_isoAng_from_energy(
533 aTrack.GetKineticEnergy(),
534 incoELM->second->find(tempLH.first / kelvin)->second);
535 anEPM_TH_E = create_E_isoAng_from_energy(
536 aTrack.GetKineticEnergy(),
537 incoELM->second->find(tempLH.second / kelvin)->second);
538 }
539 else if (tempLH.first == 0.0) {
540 // Extrapolate T0 and T1
541 anEPM_TL_E = create_E_isoAng_from_energy(
542 aTrack.GetKineticEnergy(),
543 incoELM->second->find(v_temp[0])->second);
544 anEPM_TH_E = create_E_isoAng_from_energy(
545 aTrack.GetKineticEnergy(),
546 incoELM->second->find(v_temp[1])->second);
547 tempLH.first = tempLH.second;
548 tempLH.second = v_temp[1];
549 }
550 else if (tempLH.second == 0.0) {
551 // Extrapolate Tmax-1 and Tmax
552 std::size_t nn = v_temp.size();
553 if (nn < 2) { return &theParticleChange; }
554 anEPM_TL_E = create_E_isoAng_from_energy(
555 aTrack.GetKineticEnergy(),
556 incoELM->second->find(v_temp[nn - 2])->second);
557 anEPM_TH_E = create_E_isoAng_from_energy(
558 aTrack.GetKineticEnergy(),
559 incoELM->second->find(v_temp.back())->second);
560 }
561
562 // E_isoAng for aTemp and aTrack.GetKineticEnergy()
563 G4double mu = 1.0;
564
565 // New Geant4 method - Stochastic interpolation of the final state
566 E_isoAng anEPM_T_E_sampled;
567 G4double rand_temp = G4UniformRand();
568 if (rand_temp < (aTemp - tempLH.first) / (tempLH.second - tempLH.first))
569 anEPM_T_E_sampled = std::move(anEPM_TH_E);
570 else
571 anEPM_T_E_sampled = std::move(anEPM_TL_E);
572
573 mu = getMu(&anEPM_T_E_sampled);
574
575 // Set Final State
576 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy()); // No energy change in Elastic
577 G4double phi = CLHEP::twopi * G4UniformRand();
578 G4double sint = std::sqrt(1 - mu * mu);
579 theParticleChange.SetMomentumChange(sint * std::cos(phi), sint * std::sin(phi), mu);
580 }
581 delete dp;
582
583 return &theParticleChange;
584 }
585 // Not thermal element
586 // Neutron HP will handle
587 return theHPElastic->ApplyYourself(aTrack, aNucleus,
588 true); // L. Thulliez 2021/05/04 (CEA-Saclay)
589}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
G4double GetZ() const
Definition G4Element.hh:119
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTemperature() const
const G4Element * GetElement(G4int iel) const
std::size_t GetNumberOfElements() const
G4double total(Particle const *const p1, Particle const *const p2)

◆ BuildPhysicsTable()

void G4ParticleHPThermalScattering::BuildPhysicsTable ( const G4ParticleDefinition & particle)
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 121 of file G4ParticleHPThermalScattering.cc.

122{
123 buildPhysicsTable();
124 theHPElastic->BuildPhysicsTable(particle);
125}

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4ParticleHPThermalScattering::GetFatalEnergyCheckLevels ( ) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 1182 of file G4ParticleHPThermalScattering.cc.

1183{
1184 // max energy non-conservation is mass of heavy nucleus
1185 return std::pair<G4double, G4double>(10.0 * perCent, 350.0 * CLHEP::GeV);
1186}

◆ ModelDescription()

void G4ParticleHPThermalScattering::ModelDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 1210 of file G4ParticleHPThermalScattering.cc.

1211{
1212 outFile << "High Precision model based on thermal scattering data in\n"
1213 << "evaluated nuclear data libraries for neutrons below 5eV\n"
1214 << "on specific materials\n";
1215}

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