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

#include <MCGIDI.hpp>

Inheritance diagram for MCGIDI::ProtareTNSL:

Public Member Functions

LUPI_HOST_DEVICE ProtareTNSL ()
LUPI_HOST ProtareTNSL (LUPI::StatusMessageReporting &a_smr, GIDI::ProtareTNSL const &a_protare, PoPI::Database const &a_pops, Transporting::MC &a_settings, GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash, GIDI::Styles::TemperatureInfos const &a_temperatureInfos, GIDI::ExcludeReactionsSet const &a_reactionsToExclude, std::size_t a_reactionsToExcludeOffset=0, bool a_allowFixedGrid=true)
LUPI_HOST_DEVICE ~ProtareTNSL ()
LUPI_HOST_DEVICE ProtareSingle const * protareWithElastic () const
LUPI_HOST_DEVICE ProtareSingle const * TNSL () const
LUPI_HOST_DEVICE ProtareSingle const * protareWithoutElastic () const
LUPI_HOST_DEVICE double TNSL_maximumEnergy () const
LUPI_HOST_DEVICE double TNSL_maximumTemperature () const
LUPI_HOST void setUserParticleIndex2 (int a_particleIndex, int a_userParticleIndex)
LUPI_HOST void setUserParticleIndexViaIntid2 (int a_particleIntid, int a_userParticleIndex)
LUPI_HOST_DEVICE std::size_t numberOfProtares () const
LUPI_HOST_DEVICE ProtareSingle const * protare (std::size_t a_index) const
LUPI_HOST_DEVICE ProtareSingleprotare (std::size_t a_index)
LUPI_HOST_DEVICE ProtareSingle const * protareWithReaction (std::size_t a_index) const
LUPI_HOST_DEVICE double minimumEnergy () const
LUPI_HOST_DEVICE double maximumEnergy () const
LUPI_HOST_DEVICE Vector< double > temperatures (std::size_t a_index=0) const
LUPI_HOST Vector< double > const & projectileMultiGroupBoundaries () const
LUPI_HOST Vector< double > const & projectileMultiGroupBoundariesCollapsed () const
LUPI_HOST_DEVICE std::size_t numberOfReactions () const
LUPI_HOST_DEVICE Reaction const * reaction (std::size_t a_index) const
LUPI_HOST_DEVICE std::size_t numberOfOrphanProducts () const
LUPI_HOST_DEVICE Reaction const * orphanProduct (std::size_t a_index) const
LUPI_HOST_DEVICE bool hasFission () const
LUPI_HOST_DEVICE bool hasIncoherentDoppler () const
LUPI_HOST_DEVICE int URR_index () const
LUPI_HOST_DEVICE bool hasURR_probabilityTables () const
LUPI_HOST_DEVICE double URR_domainMin () const
LUPI_HOST_DEVICE double URR_domainMax () const
LUPI_HOST_DEVICE bool reactionHasURR_probabilityTables (std::size_t a_index) const
LUPI_HOST_DEVICE double threshold (std::size_t a_index) const
LUPI_HOST_DEVICE double crossSection (URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const
LUPI_HOST_DEVICE void crossSectionVector (double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const
LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const
LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy) const
template<typename RNG>
LUPI_HOST_DEVICE std::size_t sampleReaction (Sampling::Input &a_input, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_crossSection, RNG &&a_rng) const
LUPI_HOST_DEVICE double depositionEnergy (std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double depositionMomentum (std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double productionEnergy (std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double gain (std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex) const
LUPI_HOST_DEVICE double gainViaIntid (std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid) const
LUPI_HOST_DEVICE Vector< double > const & upscatterModelAGroupVelocities () const
LUPI_HOST_DEVICE void serialize2 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE long sizeOf2 () const
Public Member Functions inherited from MCGIDI::Protare
LUPI_HOST_DEVICE Protare (ProtareType a_protareType)
LUPI_HOST Protare (ProtareType a_protareType, GIDI::Protare const &a_protare, Transporting::MC const &a_settings, PoPI::Database const &a_pops)
virtual LUPI_HOST_DEVICE ~Protare ()
LUPI_HOST_DEVICE ProtareType protareType () const
LUPI_HOST_DEVICE String const & projectileID () const
LUPI_HOST_DEVICE int projectileIntid () const
LUPI_HOST_DEVICE int projectileIndex () const
LUPI_HOST_DEVICE int projectileUserIndex () const
LUPI_HOST_DEVICE double projectileMass () const
LUPI_HOST_DEVICE double projectileExcitationEnergy () const
LUPI_HOST_DEVICE String const & targetID () const
LUPI_HOST_DEVICE int targetIntid () const
LUPI_HOST_DEVICE int targetIndex () const
LUPI_HOST_DEVICE int targetUserIndex () const
LUPI_HOST_DEVICE double targetMass () const
LUPI_HOST_DEVICE double targetExcitationEnergy () const
LUPI_HOST_DEVICE int photonIndex () const
LUPI_HOST_DEVICE int userPhotonIndex () const
LUPI_HOST_DEVICE String evaluation () const
LUPI_HOST GIDI::Frame projectileFrame () const
LUPI_HOST Vector< int > const & productIntids (bool a_transportablesOnly) const
LUPI_HOST Vector< int > const & productIndices (bool a_transportablesOnly) const
LUPI_HOST Vector< int > const & userProductIndices (bool a_transportablesOnly) const
LUPI_HOST void setUserParticleIndex (int a_particleIndex, int a_userParticleIndex)
LUPI_HOST void setUserParticleIndexViaIntid (int a_particleIntid, int a_userParticleIndex)
LUPI_HOST_DEVICE bool isTNSL_ProtareSingle () const
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfProtares () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE ProtareSingle const * protare (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE ProtareSingleprotare (std::size_t a_index) MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE ProtareSingle const * protareWithReaction (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double minimumEnergy () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double maximumEnergy () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Vector< double > temperatures (std::size_t a_index=0) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST Vector< double > const & projectileMultiGroupBoundaries () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST Vector< double > const & projectileMultiGroupBoundariesCollapsed () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfReactions () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Reaction const * reaction (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfOrphanProducts () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Reaction const * orphanProduct (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasFission () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasIncoherentDoppler () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE int URR_index () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasURR_probabilityTables () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double URR_domainMin () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double URR_domainMax () const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool reactionHasURR_probabilityTables (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double threshold (std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double crossSection (URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE void crossSectionVector (double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
template<typename RNG>
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t sampleReaction (Sampling::Input &a_input, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_crossSection, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double depositionEnergy (std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double depositionMomentum (std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double productionEnergy (std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double gain (std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double gainViaIntid (std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Vector< double > const & upscatterModelAGroupVelocities () const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize2 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serializeCommon (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE long sizeOf () const
LUPI_HOST_DEVICE long memorySize ()
LUPI_HOST_DEVICE void incrementMemorySize (long &a_totalMemory, long &a_sharedMemory)
template<typename RNG>
LUPI_HOST_DEVICE std::size_t sampleReaction (Sampling::Input &a_input, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_crossSection, RNG &&a_rng) const

Additional Inherited Members

Public Attributes inherited from MCGIDI::Protare
friend ProtareSingle
friend ProtareComposite
friend ProtareTNSL

Detailed Description

Class representing a GNDS <reactionSuite> node with only data needed for Monte Carlo transport. The data are also stored in a way that is better suited for Monte Carlo transport. For example, cross section data for each reaction are not stored with its reaction, but within the HeatedCrossSections member of the Protare.

Definition at line 1833 of file MCGIDI.hpp.

Constructor & Destructor Documentation

◆ ProtareTNSL() [1/2]

LUPI_HOST_DEVICE MCGIDI::ProtareTNSL::ProtareTNSL ( )

Default constructor used when broadcasting a Protare as needed by MPI or GPUs.

Definition at line 41 of file MCGIDI_protareTNSL.cc.

41 :
43 m_numberOfTNSLReactions( 0 ),
44 m_TNSL_maximumEnergy( 0.0 ),
45 m_TNSL_maximumTemperature( 0.0 ),
46 m_protareWithElastic( nullptr ),
47 m_TNSL( nullptr ),
48 m_protareWithoutElastic( nullptr ) {
49
50}
LUPI_HOST_DEVICE Protare(ProtareType a_protareType)

Referenced by ProtareTNSL().

◆ ProtareTNSL() [2/2]

LUPI_HOST MCGIDI::ProtareTNSL::ProtareTNSL ( LUPI::StatusMessageReporting & a_smr,
GIDI::ProtareTNSL const & a_protare,
PoPI::Database const & a_pops,
Transporting::MC & a_settings,
GIDI::Transporting::Particles const & a_particles,
DomainHash const & a_domainHash,
GIDI::Styles::TemperatureInfos const & a_temperatureInfos,
GIDI::ExcludeReactionsSet const & a_reactionsToExclude,
std::size_t a_reactionsToExcludeOffset = 0,
bool a_allowFixedGrid = true )

◆ ~ProtareTNSL()

LUPI_HOST_DEVICE MCGIDI::ProtareTNSL::~ProtareTNSL ( )

Definition at line 100 of file MCGIDI_protareTNSL.cc.

100 {
101
102 delete m_protareWithElastic;
103 delete m_TNSL;
104 delete m_protareWithoutElastic;
105}

Member Function Documentation

◆ crossSection()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::crossSection ( URR_protareInfos const & a_URR_protareInfos,
std::size_t a_hashIndex,
double a_temperature,
double a_energy,
bool a_sampling = false ) const

Returns the total cross section.

Parameters
a_URR_protareInfos[in] URR information.
a_hashIndex[in] The cross section hash index.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
a_sampling[in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold is present the cross section is augmented.
Returns
The total cross section.

Definition at line 256 of file MCGIDI_protareTNSL.cc.

256 {
257
258 double crossSection1 = 0.0;
259
260 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
261 crossSection1 = m_TNSL->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling ) +
262 m_protareWithoutElastic->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling ); }
263 else {
264 crossSection1 = m_protareWithElastic->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
265 }
266
267 return( crossSection1 );
268}

◆ crossSectionVector()

LUPI_HOST_DEVICE void MCGIDI::ProtareTNSL::crossSectionVector ( double a_temperature,
double a_userFactor,
std::size_t a_numberAllocated,
double * a_crossSectionVector ) const

Adds the energy dependent, total cross section corresponding to the temperature a_temperature multiplied by a_userFactor to a_crossSectionVector.

Parameters
a_temperature[in] Specifies the temperature of the material.
a_userFactor[in] User factor which all cross sections are multiplied by.
a_numberAllocated[in] The length of memory allocated for a_crossSectionVector.
a_crossSectionVector[in/out] The energy dependent, total cross section to add cross section data to.

Definition at line 279 of file MCGIDI_protareTNSL.cc.

280 {
281
282 if( a_temperature <= m_TNSL_maximumTemperature ) {
283 m_TNSL->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
284 m_protareWithoutElastic->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector ); }
285 else {
286 m_protareWithElastic->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
287 }
288}

◆ depositionEnergy()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::depositionEnergy ( std::size_t a_hashIndex,
double a_temperature,
double a_energy ) const

Returns the total deposition energy.

Parameters
a_hashIndex[in] The cross section hash index.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
Returns
The total deposition energy.

Definition at line 361 of file MCGIDI_protareTNSL.cc.

361 {
362
363 double deposition_energy = 0.0;
364
365 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
366 deposition_energy = m_TNSL->depositionEnergy( a_hashIndex, a_temperature, a_energy ) +
367 m_protareWithoutElastic->depositionEnergy( a_hashIndex, a_temperature, a_energy ); }
368 else {
369 deposition_energy = m_protareWithElastic->depositionEnergy( a_hashIndex, a_temperature, a_energy );
370 }
371
372 return( deposition_energy );
373}

◆ depositionMomentum()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::depositionMomentum ( std::size_t a_hashIndex,
double a_temperature,
double a_energy ) const

Returns the total deposition momentum.

Parameters
a_hashIndex[in] The cross section hash index.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
Returns
The total deposition momentum.

Definition at line 385 of file MCGIDI_protareTNSL.cc.

385 {
386
387 double deposition_momentum = 0.0;
388
389 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
390 deposition_momentum = m_TNSL->depositionMomentum( a_hashIndex, a_temperature, a_energy ) +
391 m_protareWithoutElastic->depositionMomentum( a_hashIndex, a_temperature, a_energy ); }
392 else {
393 deposition_momentum = m_protareWithElastic->depositionMomentum( a_hashIndex, a_temperature, a_energy );
394 }
395
396 return( deposition_momentum );
397}

◆ gain()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::gain ( std::size_t a_hashIndex,
double a_temperature,
double a_energy,
int a_particleIndex ) const

Returns the gain for particle with index a_particleIndex.

Parameters
a_hashIndex[in] The cross section hash index.
a_temperature[in] The temperature of the target.
a_energy[in] The projectile energy.
a_particleIndex[in] The index of the particle whose gain is to be returned.
Returns
[in] A vector of the length of the number of multi-group groups.

Definition at line 434 of file MCGIDI_protareTNSL.cc.

434 {
435
436 double gain1 = 0.0;
437
438 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
439 gain1 = m_TNSL->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex ) +
440 m_protareWithoutElastic->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex ); }
441 else {
442 gain1 = m_protareWithElastic->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
443 }
444
445 return( gain1 );
446}

◆ gainViaIntid()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::gainViaIntid ( std::size_t a_hashIndex,
double a_temperature,
double a_energy,
int a_particleIntid ) const

Returns the gain for particle with intid a_particleIntid.

Parameters
a_hashIndex[in] The cross section hash index.
a_temperature[in] The temperature of the target.
a_energy[in] The projectile energy.
a_particleIntid[in] The intid of the particle whose gain is to be returned.
Returns
[in] A vector of the length of the number of multi-group groups.

Definition at line 459 of file MCGIDI_protareTNSL.cc.

459 {
460
461 double gain1 = 0.0;
462
463 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
464 gain1 = m_TNSL->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid ) +
465 m_protareWithoutElastic->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid ); }
466 else {
467 gain1 = m_protareWithElastic->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
468 }
469
470 return( gain1 );
471}

◆ hasFission()

LUPI_HOST_DEVICE bool MCGIDI::ProtareTNSL::hasFission ( ) const
inline

Definition at line 1882 of file MCGIDI.hpp.

1882{ return( m_protareWithElastic->hasFission( ) ); } /* Returns the normal ProtareSingle's hasFission value. */

◆ hasIncoherentDoppler()

LUPI_HOST_DEVICE bool MCGIDI::ProtareTNSL::hasIncoherentDoppler ( ) const
inline

Definition at line 1883 of file MCGIDI.hpp.

1883{ return( false ); } /* Always returns false as this is a neutron as projectile and not a photon. */

◆ hasURR_probabilityTables()

LUPI_HOST_DEVICE bool MCGIDI::ProtareTNSL::hasURR_probabilityTables ( ) const
inline

Definition at line 1886 of file MCGIDI.hpp.

1886{ return( m_protareWithElastic->hasURR_probabilityTables( ) ); }

◆ maximumEnergy()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::maximumEnergy ( ) const
inline

Returns the maximum cross section domain.

Definition at line 1867 of file MCGIDI.hpp.

◆ minimumEnergy()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::minimumEnergy ( ) const
inline

Returns the minimum cross section domain.

Definition at line 1866 of file MCGIDI.hpp.

◆ numberOfOrphanProducts()

LUPI_HOST_DEVICE std::size_t MCGIDI::ProtareTNSL::numberOfOrphanProducts ( ) const
inline

Returns the number of orphan products in the normal ProtareSingle.

Definition at line 1877 of file MCGIDI.hpp.

◆ numberOfProtares()

LUPI_HOST_DEVICE std::size_t MCGIDI::ProtareTNSL::numberOfProtares ( ) const
inline

Always Returns 2.

Definition at line 1861 of file MCGIDI.hpp.

◆ numberOfReactions()

LUPI_HOST_DEVICE std::size_t MCGIDI::ProtareTNSL::numberOfReactions ( ) const
inline

Definition at line 1875 of file MCGIDI.hpp.

1875{ return( m_TNSL->numberOfReactions( ) + m_protareWithElastic->numberOfReactions( ) ); }

◆ orphanProduct()

LUPI_HOST_DEVICE Reaction const * MCGIDI::ProtareTNSL::orphanProduct ( std::size_t a_index) const
inline

Returns the (a_index - 1 )^th orphan product in the normal ProtareSingle.

Definition at line 1879 of file MCGIDI.hpp.

◆ productionEnergy()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::productionEnergy ( std::size_t a_hashIndex,
double a_temperature,
double a_energy ) const

Returns the total production energy.

Parameters
a_hashIndex[in] The cross section hash index.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
Returns
The total production energy.

Definition at line 409 of file MCGIDI_protareTNSL.cc.

409 {
410
411 double production_energy = 0.0;
412
413 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
414 production_energy = m_TNSL->productionEnergy( a_hashIndex, a_temperature, a_energy ) +
415 m_protareWithoutElastic->productionEnergy( a_hashIndex, a_temperature, a_energy ); }
416 else {
417 production_energy = m_protareWithElastic->productionEnergy( a_hashIndex, a_temperature, a_energy );
418 }
419
420 return( production_energy );
421}

◆ projectileMultiGroupBoundaries()

LUPI_HOST Vector< double > const & MCGIDI::ProtareTNSL::projectileMultiGroupBoundaries ( ) const
inline

Returns the value of the m_projectileMultiGroupBoundaries member.

Definition at line 1870 of file MCGIDI.hpp.

◆ projectileMultiGroupBoundariesCollapsed()

LUPI_HOST Vector< double > const & MCGIDI::ProtareTNSL::projectileMultiGroupBoundariesCollapsed ( ) const
inline

Returns the value of the m_projectileMultiGroupBoundariesCollapsed member.

Definition at line 1872 of file MCGIDI.hpp.

◆ protare() [1/2]

LUPI_HOST_DEVICE ProtareSingle * MCGIDI::ProtareTNSL::protare ( std::size_t a_index)

Returns the pointer representing the (a_index - 1)th ProtareSingle.

Parameters
a_index[in] Index of the ProtareSingle to return. Can only be 0 or 1.
Returns
Pointer to the requested protare or nullptr if invalid a_index..

Definition at line 159 of file MCGIDI_protareTNSL.cc.

159 {
160
161 if( a_index == 0 ) return( m_protareWithElastic );
162 if( a_index == 1 ) return( m_TNSL );
163 return( nullptr );
164}

◆ protare() [2/2]

LUPI_HOST_DEVICE ProtareSingle const * MCGIDI::ProtareTNSL::protare ( std::size_t a_index) const

Returns the pointer representing the (a_index - 1)th ProtareSingle.

Parameters
a_index[in] Index of the ProtareSingle to return. Can only be 0 or 1.
Returns
Pointer to the requested protare or nullptr if invalid a_index..

Definition at line 144 of file MCGIDI_protareTNSL.cc.

144 {
145
146 if( a_index == 0 ) return( m_protareWithElastic );
147 if( a_index == 1 ) return( m_TNSL );
148 return( nullptr );
149}

◆ protareWithElastic()

LUPI_HOST_DEVICE ProtareSingle const * MCGIDI::ProtareTNSL::protareWithElastic ( ) const
inline

Returns the m_protareWithElastic member.

Definition at line 1850 of file MCGIDI.hpp.

◆ protareWithoutElastic()

LUPI_HOST_DEVICE ProtareSingle const * MCGIDI::ProtareTNSL::protareWithoutElastic ( ) const
inline

Returns the m_protareWithoutElastic member.

Definition at line 1852 of file MCGIDI.hpp.

◆ protareWithReaction()

LUPI_HOST_DEVICE ProtareSingle const * MCGIDI::ProtareTNSL::protareWithReaction ( std::size_t a_index) const

Returns the pointer to the ProtareSingle that contains the (a_index - 1)th reaction.

Parameters
a_index[in] Index of the reaction.
Returns
Pointer to the requested protare or nullptr if invalid a_index..

Definition at line 174 of file MCGIDI_protareTNSL.cc.

174 {
175
176 if( a_index < m_numberOfTNSLReactions ) return( m_TNSL );
177 return( m_protareWithElastic->protareWithReaction( a_index - m_numberOfTNSLReactions ) );
178}

◆ reaction()

LUPI_HOST_DEVICE Reaction const * MCGIDI::ProtareTNSL::reaction ( std::size_t a_index) const

Returns the reaction at index a_index. If a_index is negative, the reaction of the TNSL protare at index -*a_index* is returned; otherwise, the reaction from the regular protare at index a_index is returned.

Parameters
a_index[in] The index of the reaction to return.
Returns
The reaction at index a_index.

Definition at line 208 of file MCGIDI_protareTNSL.cc.

208 {
209
210 if( a_index < m_numberOfTNSLReactions ) return( m_TNSL->reaction( a_index ) );
211 return( m_protareWithElastic->reaction( a_index - m_numberOfTNSLReactions ) );
212}

◆ reactionCrossSection() [1/2]

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::reactionCrossSection ( std::size_t a_reactionIndex,
URR_protareInfos const & a_URR_protareInfos,
double a_temperature,
double a_energy ) const

Returns the cross section for reaction at index a_reactionIndex, for target at temperature a_temperature and projectile of energy a_energy.

Parameters
a_URR_protareInfos[in] URR information.
a_reactionIndex[in] The index of the reaction.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
Returns
The total cross section.

Definition at line 333 of file MCGIDI_protareTNSL.cc.

333 {
334
335 std::size_t index = a_reactionIndex - m_numberOfTNSLReactions;
336 double crossSection1 = 0.0;
337
338 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
339 if( a_reactionIndex < m_numberOfTNSLReactions ) {
340 crossSection1 = m_TNSL->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_temperature, a_energy ); }
341 else {
342 if( a_reactionIndex > m_numberOfTNSLReactions ) crossSection1 = m_protareWithElastic->reactionCrossSection( index, a_URR_protareInfos, a_temperature, a_energy );
343 } }
344 else {
345 if( a_reactionIndex >= m_numberOfTNSLReactions ) crossSection1 = m_protareWithElastic->reactionCrossSection( index, a_URR_protareInfos, a_temperature, a_energy );
346 }
347
348 return( crossSection1 );
349}

◆ reactionCrossSection() [2/2]

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::reactionCrossSection ( std::size_t a_reactionIndex,
URR_protareInfos const & a_URR_protareInfos,
std::size_t a_hashIndex,
double a_temperature,
double a_energy,
bool a_sampling = false ) const

Returns the cross section for reaction at index a_reactionIndex, for target at temperature a_temperature and projectile of energy a_energy.

Parameters
a_URR_protareInfos[in] URR information.
a_reactionIndex[in] The index of the reaction.
a_hashIndex[in] The cross section hash index.
a_temperature[in] The target temperature.
a_energy[in] The projectile energy.
a_sampling[in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold is present the cross section is augmented.
Returns
The total cross section.

Definition at line 304 of file MCGIDI_protareTNSL.cc.

304 {
305
306 std::size_t index = a_reactionIndex - m_numberOfTNSLReactions;
307 double crossSection1 = 0.0;
308
309 if( ( a_energy < m_TNSL_maximumEnergy ) && ( a_temperature <= m_TNSL_maximumTemperature ) ) {
310 if( a_reactionIndex < m_numberOfTNSLReactions ) {
311 crossSection1 = m_TNSL->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling ); }
312 else {
313 if( a_reactionIndex > m_numberOfTNSLReactions ) crossSection1 = m_protareWithElastic->reactionCrossSection( index, a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
314 } }
315 else {
316 if( a_reactionIndex >= m_numberOfTNSLReactions ) crossSection1 = m_protareWithElastic->reactionCrossSection( index, a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
317 }
318
319 return( crossSection1 );
320}

◆ reactionHasURR_probabilityTables()

LUPI_HOST_DEVICE bool MCGIDI::ProtareTNSL::reactionHasURR_probabilityTables ( std::size_t a_index) const

Returns true if the reaction at index a_index has URR robability tables and false otherwise.

Parameters
a_index[in] The index of the reaction.
Returns
true if the reaction has URR robability tables and false otherwise.

Definition at line 222 of file MCGIDI_protareTNSL.cc.

222 {
223
224 if( a_index < m_numberOfTNSLReactions ) return( false );
225 return( m_protareWithElastic->reactionHasURR_probabilityTables( a_index - m_numberOfTNSLReactions ) );
226}

◆ sampleReaction()

template<typename RNG>
LUPI_HOST_DEVICE std::size_t MCGIDI::ProtareTNSL::sampleReaction ( Sampling::Input & a_input,
URR_protareInfos const & a_URR_protareInfos,
std::size_t a_hashIndex,
double a_crossSection,
RNG && a_rng ) const
inline

Returns the total cross section.

Parameters
a_input[in] Sample options requested by user.
a_URR_protareInfos[in] URR information.
a_hashIndex[in] Specifies the continuous energy hash index or multi-group index.
a_crossSection[in] The total cross section for the protare at a_input.temperature() and a_input.energy().
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
Returns
The index of the sampled reaction.

Definition at line 2929 of file MCGIDI_headerSource.hpp.

2930 {
2931
2932 std::size_t reactionIndex = 0;
2933
2934 if( ( a_input.energy( ) < m_TNSL_maximumEnergy ) && ( a_input.temperature( ) <= m_TNSL_maximumTemperature ) ) {
2935 double TNSL_crossSection = m_TNSL->crossSection( a_URR_protareInfos, a_hashIndex, a_input.temperature( ), a_input.energy( ), true );
2936
2937 if( TNSL_crossSection > a_rng( ) * a_crossSection ) {
2938 reactionIndex = m_TNSL->sampleReaction( a_input, a_URR_protareInfos, a_hashIndex, TNSL_crossSection, a_rng ); }
2939 else {
2940 reactionIndex = m_protareWithoutElastic->sampleReaction( a_input, a_URR_protareInfos, a_hashIndex,
2941 a_crossSection - TNSL_crossSection, a_rng );
2942 if( reactionIndex != MCGIDI_nullReaction ) reactionIndex += m_numberOfTNSLReactions + 1;
2943 } }
2944 else {
2945 reactionIndex = m_protareWithElastic->sampleReaction( a_input, a_URR_protareInfos, a_hashIndex, a_crossSection, a_rng );
2946 if( reactionIndex != MCGIDI_nullReaction ) reactionIndex += m_numberOfTNSLReactions;
2947 }
2948
2949 return( reactionIndex );
2950}
#define MCGIDI_nullReaction
Definition MCGIDI.hpp:64

◆ serialize2()

LUPI_HOST_DEVICE void MCGIDI::ProtareTNSL::serialize2 ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )

This method serializes this for broadcasting as needed for MPI and GPUs. The method can count the number of required bytes, pack this or unpack this depending on a_mode.

Parameters
a_buffer[in] The buffer to read or write data to depending on a_mode.
a_mode[in] Specifies the action of this method.

Definition at line 481 of file MCGIDI_protareTNSL.cc.

481 {
482
483 int numberOfTNSLReactions = static_cast<int>( m_numberOfTNSLReactions );
484 DATA_MEMBER_INT( numberOfTNSLReactions, a_buffer, a_mode );
485 m_numberOfTNSLReactions = static_cast<std::size_t>( numberOfTNSLReactions );
486
487 DATA_MEMBER_DOUBLE( m_TNSL_maximumEnergy, a_buffer, a_mode );
488 DATA_MEMBER_DOUBLE( m_TNSL_maximumTemperature, a_buffer, a_mode );
489
490 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
491 if( a_buffer.m_placement != nullptr ) {
492 m_protareWithElastic = new(a_buffer.m_placement) ProtareSingle;
493 a_buffer.incrementPlacement( sizeof( ProtareSingle ) );
494 m_TNSL = new(a_buffer.m_placement) ProtareSingle;
495 a_buffer.incrementPlacement( sizeof( ProtareSingle ) );
496 m_protareWithoutElastic = new(a_buffer.m_placement) ProtareSingle;
497 a_buffer.incrementPlacement( sizeof( ProtareSingle ) ); }
498 else {
499 m_protareWithElastic = new ProtareSingle( );
500 m_TNSL = new ProtareSingle( );
501 m_protareWithoutElastic = new ProtareSingle( );
502 }
503 }
504 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
505 a_buffer.incrementPlacement( sizeof( ProtareSingle ) );
506 a_buffer.incrementPlacement( sizeof( ProtareSingle ) );
507 a_buffer.incrementPlacement( sizeof( ProtareSingle ) );
508 }
509 m_protareWithElastic->serializeCommon( a_buffer, a_mode );
510 m_protareWithElastic->serialize2( a_buffer, a_mode );
511 m_TNSL->serializeCommon( a_buffer, a_mode );
512 m_TNSL->serialize2( a_buffer, a_mode );
513 m_protareWithoutElastic->serializeCommon( a_buffer, a_mode );
514 m_protareWithoutElastic->serialize2( a_buffer, a_mode );
515}
#define DATA_MEMBER_DOUBLE(member, buf, mode)
#define DATA_MEMBER_INT( member, buf, mode)
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
friend ProtareSingle
Definition MCGIDI.hpp:1599

◆ setUserParticleIndex2()

LUPI_HOST void MCGIDI::ProtareTNSL::setUserParticleIndex2 ( int a_particleIndex,
int a_userParticleIndex )

Updates the m_userParticleIndex to a_userParticleIndex for all particles with PoPs index a_particleIndex.

Parameters
a_particleIndex[in] The PoPs index of the particle whose user index is to be set.
a_userParticleIndex[in] The particle index specified by the user.

Definition at line 114 of file MCGIDI_protareTNSL.cc.

114 {
115
116 m_protareWithElastic->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
117 m_TNSL->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
118 m_protareWithoutElastic->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
119}

◆ setUserParticleIndexViaIntid2()

LUPI_HOST void MCGIDI::ProtareTNSL::setUserParticleIndexViaIntid2 ( int a_particleIntid,
int a_userParticleIndex )

Updates the m_userParticleIndex to a_userParticleIndex for all particles with PoPs intid a_particleIntid.

Parameters
a_particleIntid[in] The PoPs intid of the particle whose user index is to be set.
a_userParticleIndex[in] The particle index specified by the user.

Definition at line 128 of file MCGIDI_protareTNSL.cc.

128 {
129
130 m_protareWithElastic->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
131 m_TNSL->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
132 m_protareWithoutElastic->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
133}

◆ sizeOf2()

LUPI_HOST_DEVICE long MCGIDI::ProtareTNSL::sizeOf2 ( ) const
inline

Definition at line 1914 of file MCGIDI.hpp.

1914{ return sizeof(*this); }

◆ temperatures()

LUPI_HOST_DEVICE Vector< double > MCGIDI::ProtareTNSL::temperatures ( std::size_t a_index = 0) const

Returns the list of temperatures for the requested ProtareSingle.

Parameters
a_index[in] Index of the reqested ProtareSingle.
Returns
Vector of doubles.

Definition at line 188 of file MCGIDI_protareTNSL.cc.

188 {
189
190 if( a_index == 0 ) return( m_protareWithElastic->temperatures( 0 ) );
191 if( a_index == 1 ) return( m_TNSL->temperatures( 0 ) );
192
193 LUPI_THROW( "ProtareSingle::temperatures: a_index not 0 or 1." );
194
195 Vector<double> temps; // Only to stop compilers from complaining.
196 return( temps );
197}
#define LUPI_THROW(arg)

◆ threshold()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::threshold ( std::size_t a_index) const

Returns the threshold for the reaction at index a_index. If a_index is negative, it is set to 0 before the threshold in the regular protare is returned.

Parameters
a_index[in] The index of the reaction.
Returns
The threshold for reaction at index a_index.

Definition at line 237 of file MCGIDI_protareTNSL.cc.

237 {
238
239 if( a_index < m_numberOfTNSLReactions ) return( m_TNSL->threshold( a_index ) );
240 return( m_protareWithElastic->threshold( a_index - m_numberOfTNSLReactions ) );
241}

◆ TNSL()

LUPI_HOST_DEVICE ProtareSingle const * MCGIDI::ProtareTNSL::TNSL ( ) const
inline

Returns the m_TNSL member.

Definition at line 1851 of file MCGIDI.hpp.

Referenced by ProtareTNSL().

◆ TNSL_maximumEnergy()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::TNSL_maximumEnergy ( ) const
inline

Definition at line 1854 of file MCGIDI.hpp.

1854{ return( m_TNSL_maximumEnergy ); }

◆ TNSL_maximumTemperature()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::TNSL_maximumTemperature ( ) const
inline

Definition at line 1855 of file MCGIDI.hpp.

1855{ return( m_TNSL_maximumTemperature ); }

◆ upscatterModelAGroupVelocities()

LUPI_HOST_DEVICE Vector< double > const & MCGIDI::ProtareTNSL::upscatterModelAGroupVelocities ( ) const
inline

Returns a reference to the m_upscatterModelAGroupVelocities member.

Definition at line 1910 of file MCGIDI.hpp.

◆ URR_domainMax()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::URR_domainMax ( ) const
inline

Definition at line 1888 of file MCGIDI.hpp.

1888{ return( m_protareWithElastic->URR_domainMax( ) ); }

◆ URR_domainMin()

LUPI_HOST_DEVICE double MCGIDI::ProtareTNSL::URR_domainMin ( ) const
inline

Definition at line 1887 of file MCGIDI.hpp.

1887{ return( m_protareWithElastic->URR_domainMin( ) ); }

◆ URR_index()

LUPI_HOST_DEVICE int MCGIDI::ProtareTNSL::URR_index ( ) const
inline

Definition at line 1885 of file MCGIDI.hpp.

1885{ return( -1 ); }

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