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

#include <MCGIDI.hpp>

Public Member Functions

LUPI_HOST_DEVICE Reaction ()
LUPI_HOST Reaction (GIDI::Reaction const &a_reaction, SetupInfo &a_setupInfo, Transporting::MC const &a_settings, GIDI::Transporting::Particles const &a_particles, GIDI::Styles::TemperatureInfos const &a_temperatureInfos)
LUPI_HOST_DEVICE ~Reaction ()
LUPI_HOST_DEVICE void updateProtareSingleInfo (ProtareSingle *a_protareSingle, std::size_t a_reactionIndex)
LUPI_HOST_DEVICE ProtareSingle const * protareSingle () const
LUPI_HOST_DEVICE std::size_t reactionIndex () const
LUPI_HOST_DEVICE std::size_t GIDI_reactionIndex () const
LUPI_HOST_DEVICE String const & label () const
LUPI_HOST_DEVICE int ENDF_MT () const
LUPI_HOST_DEVICE int ENDL_C () const
LUPI_HOST_DEVICE int ENDL_S () const
LUPI_HOST_DEVICE int initialStateIndex () const
LUPI_HOST_DEVICE double finalQ (double a_energy) const
LUPI_HOST_DEVICE bool hasFission () const
LUPI_HOST_DEVICE double projectileMass () const
LUPI_HOST_DEVICE double targetMass () const
LUPI_HOST_DEVICE double crossSectionThreshold () const
LUPI_HOST_DEVICE double twoBodyThreshold () const
LUPI_HOST_DEVICE double crossSection (URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double crossSection (URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy) const
LUPI_HOST GIDI::Functions::XYs1d crossSectionAsGIDI_XYs1d (double a_temperature) const
LUPI_HOST_DEVICE Vector< int > const & productIntids () const
LUPI_HOST_DEVICE Vector< int > const & productIndices () const
LUPI_HOST_DEVICE Vector< int > const & userProductIndices () const
LUPI_HOST_DEVICE std::size_t numberOfProducts () const
LUPI_HOST_DEVICE Product const * product (std::size_t a_index) const
LUPI_HOST_DEVICE int productMultiplicity (int a_index) const
LUPI_HOST_DEVICE int productMultiplicityViaIntid (int a_intid) const
LUPI_HOST_DEVICE int productMultiplicities (int a_index) const
LUPI_HOST_DEVICE double productAverageMultiplicity (int a_index, double a_projectileEnergy) const
LUPI_HOST_DEVICE double productAverageMultiplicityViaIntid (int a_intid, double a_projectileEnergy) const
LUPI_HOST_DEVICE Vector< int > const & productIntidsTransportable () const
LUPI_HOST_DEVICE Vector< int > const & productIndicesTransportable () const
LUPI_HOST_DEVICE Vector< int > const & userProductIndicesTransportable () const
LUPI_HOST_DEVICE Vector< std::size_t > associatedOrphanProductIndices () const
LUPI_HOST void addOrphanProductToProductList (std::vector< Product * > &a_associatedOrphanProducts) const
LUPI_HOST_DEVICE void addOrphanProductToProductList (Vector< Product * > &a_associatedOrphanProducts) const
LUPI_HOST_DEVICE void addOrphanProductToProductList (Vector< Reaction * > &a_orphanProducts)
LUPI_HOST void setOrphanProductData (std::vector< std::size_t > const &a_associatedOrphanProductIndcies, std::vector< Product * > const &a_associatedOrphanProducts)
LUPI_HOST void setUserParticleIndex (int a_particleIndex, int a_userParticleIndex)
LUPI_HOST void setUserParticleIndexViaIntid (int a_particleIntid, int a_userParticleIndex)
LUPI_HOST void setModelDBRC_data (Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void sampleProducts (Protare const *a_protare, Sampling::Input &a_input, RNG &&a_rng, PUSHBACK &&a_push_back, Sampling::ProductHandler &a_products, bool a_checkOrphanProducts=true) const
template<typename RNG>
LUPI_HOST_DEVICE double angleBiasing (int a_pid, double a_temperature, double a_energy_in, double a_mu_lab, double &a_energy_out, RNG &&a_rng, double *a_cumulative_weight=nullptr, bool a_checkOrphanProducts=true) const
template<typename RNG>
LUPI_HOST_DEVICE double angleBiasingViaIntid (int a_intid, double a_temperature, double a_energy_in, double a_mu_lab, double &a_energy_out, RNG &&a_rng, double *a_cumulative_weight=nullptr, bool a_checkOrphanProducts=true) const
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

Static Public Member Functions

template<typename RNG, typename PUSHBACK>
static LUPI_HOST_DEVICE void sampleNullProducts (Protare const &a_protare, double a_projectileEnergy, Sampling::Input &a_input, RNG &&a_rng, PUSHBACK &&a_push_back, Sampling::ProductHandler &a_products)

Detailed Description

Class representing a GNDS <reaction> node with only data needed for Monte Carlo transport.

Definition at line 1335 of file MCGIDI.hpp.

Constructor & Destructor Documentation

◆ Reaction() [1/2]

LUPI_HOST_DEVICE MCGIDI::Reaction::Reaction ( )

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

Definition at line 22 of file MCGIDI_reaction.cc.

22 :
23 m_protareSingle( nullptr ),
24 m_reactionIndex( MCGIDI_nullReaction ),
25 m_GIDI_reactionIndex( MCGIDI_nullReaction ),
26 m_label( ),
27 m_ENDF_MT( 0 ),
28 m_ENDL_C( 0 ),
29 m_ENDL_S( 0 ),
30 m_initialStateIndex( -1 ),
31 m_neutronIndex( -1 ),
32 m_hasFission( false ),
33 m_projectileMass( 0.0 ),
34 m_targetMass( 0.0 ),
35 m_crossSectionThreshold( 0.0 ),
36 m_twoBodyThreshold( 0.0 ),
37 m_hasFinalStatePhotons( false ),
38 m_fissionResiduaIntid( -1 ),
39 m_fissionResiduaIndex( -1 ),
40 m_fissionResiduaUserIndex( -1 ),
42 m_fissionResidualMass( 0.0 ),
43 m_totalDelayedNeutronMultiplicity( nullptr ),
44#ifdef MCGIDI_USE_OUTPUT_CHANNEL
45 m_outputChannel( nullptr ),
46#endif
47
48// GRIN extras:
49 m_GRIN_specialSampleProducts( false ),
50 m_GRIN_inelasticThreshold( 0.0 ),
51 m_GRIN_maximumCaptureIncidentEnergy( 0.0 ),
52 m_GRIN_inelastic( nullptr ),
53 m_GRIN_capture( nullptr ) {
54}
#define MCGIDI_nullReaction
Definition MCGIDI.hpp:64

Referenced by addOrphanProductToProductList().

◆ Reaction() [2/2]

LUPI_HOST MCGIDI::Reaction::Reaction ( GIDI::Reaction const & a_reaction,
SetupInfo & a_setupInfo,
Transporting::MC const & a_settings,
GIDI::Transporting::Particles const & a_particles,
GIDI::Styles::TemperatureInfos const & a_temperatureInfos )

◆ ~Reaction()

LUPI_HOST_DEVICE MCGIDI::Reaction::~Reaction ( )

Definition at line 178 of file MCGIDI_reaction.cc.

178 {
179
180#ifdef MCGIDI_USE_OUTPUT_CHANNEL
181 delete m_outputChannel;
182#else
183 delete m_totalDelayedNeutronMultiplicity;
184 for( auto iter = m_products.begin( ); iter != m_products.end( ); ++iter ) delete *iter;
185 for( auto iter = m_delayedNeutrons.begin( ); iter != m_delayedNeutrons.end( ); ++iter ) delete *iter;
186 for( auto iter = m_Qs.begin( ); iter != m_Qs.end( ); ++iter ) delete *iter;
187#endif
188}

Member Function Documentation

◆ addOrphanProductToProductList() [1/3]

LUPI_HOST void MCGIDI::Reaction::addOrphanProductToProductList ( std::vector< Product * > & a_associatedOrphanProducts) const

Adds the associated orphan products of an orphan product reaction to a_associatedOrphanProducts.

Parameters
a_associatedOrphanProducts[in] A list where the associated orphan products are added to.

Definition at line 460 of file MCGIDI_reaction.cc.

460 {
461
462#ifdef MCGIDI_USE_OUTPUT_CHANNEL
463 m_outputChannel->addOrphanProductToProductList( a_associatedOrphanProducts );
464#else
465 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
466 a_associatedOrphanProducts.push_back( *productIter );
467 }
468#endif
469}

Referenced by addOrphanProductToProductList().

◆ addOrphanProductToProductList() [2/3]

LUPI_HOST_DEVICE void MCGIDI::Reaction::addOrphanProductToProductList ( Vector< Product * > & a_associatedOrphanProducts) const

Adds the associated orphan products of an orphan product reaction to a_associatedOrphanProducts.

Parameters
a_associatedOrphanProducts[in] A list where the associated orphan products are added to.

Definition at line 477 of file MCGIDI_reaction.cc.

477 {
478
479#ifdef MCGIDI_USE_OUTPUT_CHANNEL
480 m_outputChannel->addOrphanProductToProductList( a_associatedOrphanProducts );
481#else
482 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
483 a_associatedOrphanProducts.push_back( *productIter );
484 }
485#endif
486}

◆ addOrphanProductToProductList() [3/3]

LUPI_HOST_DEVICE void MCGIDI::Reaction::addOrphanProductToProductList ( Vector< Reaction * > & a_orphanProducts)

Adds the associated orphan products of an orphan product to a_associatedOrphanProducts.

Parameters
a_associatedOrphanProducts[in] A list where the associated orphan products are added to.

Definition at line 494 of file MCGIDI_reaction.cc.

494 {
495
496 for( auto associatedOrphanProductIndex = m_associatedOrphanProductIndices.begin( );
497 associatedOrphanProductIndex != m_associatedOrphanProductIndices.end( ); ++associatedOrphanProductIndex ) {
498 Reaction *orphanProduct = a_orphanProducts[*associatedOrphanProductIndex];
499 orphanProduct->addOrphanProductToProductList( m_associatedOrphanProducts );
500 }
501}
LUPI_HOST_DEVICE Reaction()

◆ angleBiasing()

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Reaction::angleBiasing ( int a_index,
double a_temperature,
double a_energy_in,
double a_mu_lab,
double & a_energy_out,
RNG && a_rng,
double * a_cumulative_weight = nullptr,
bool a_checkOrphanProducts = true ) const
inline

Returns the weight for a project with energy a_energy_in to cause this reaction to emitted a particle of index a_index at angle a_mu_lab as seen in the lab frame. If a particle is emitted, a_energy_out is its sampled outgoing energy.

Parameters
a_index[in] The index of the particle to emit.
a_temperature[in] Specifies the temperature of the material.
a_energy_in[in] The energy of the incident particle.
a_mu_lab[in] The desired mu in the lab frame for the emitted particle.
a_energy_out[in] The energy of the emitted outgoing particle.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_cumulative_weight[in] The cumulative multiplicity.
a_checkOrphanProducts[in] If true, associated orphan products are also sampled.
Returns
The weight that the particle is emitted into mu a_mu_lab.

Definition at line 3274 of file MCGIDI_headerSource.hpp.

3275 {
3276
3277 double cumulative_weight1 = 0.0;
3278 if( a_cumulative_weight == nullptr ) a_cumulative_weight = &cumulative_weight1;
3279 double weight1 = 0.0;
3280
3281#ifdef MCGIDI_USE_OUTPUT_CHANNEL
3282 m_outputChannel->angleBiasing( this, a_index, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3283#else
3284
3285 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
3286 (*productIter)->angleBiasing( this, a_index, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3287 }
3288
3289 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_index == m_neutronIndex ) ) {
3290 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
3291 DelayedNeutron const *delayedNeutron1 = m_delayedNeutrons[i1];
3292 Product const &product = delayedNeutron1->product( );
3293
3294 product.angleBiasing( this, a_index, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3295 }
3296 }
3297#endif
3298
3299 if( a_checkOrphanProducts ) {
3300 for( auto productIter = m_associatedOrphanProducts.begin( ); productIter != m_associatedOrphanProducts.end( ); ++productIter ) {
3301 (*productIter)->angleBiasing( this, a_index, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng,
3302 *a_cumulative_weight );
3303 }
3304 }
3305
3306 return( weight1 );
3307}
LUPI_HOST_DEVICE Product const * product(std::size_t a_index) const
Definition MCGIDI.hpp:1422

◆ angleBiasingViaIntid()

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Reaction::angleBiasingViaIntid ( int a_intid,
double a_temperature,
double a_energy_in,
double a_mu_lab,
double & a_energy_out,
RNG && a_rng,
double * a_cumulative_weight = nullptr,
bool a_checkOrphanProducts = true ) const
inline

Returns the weight for a project with energy a_energy_in to cause this reaction to emitted a particle of intid a_intid at angle a_mu_lab as seen in the lab frame. If a particle is emitted, a_energy_out is its sampled outgoing energy.

Parameters
a_intid[in] The intid of the particle to emit.
a_temperature[in] Specifies the temperature of the material.
a_energy_in[in] The energy of the incident particle.
a_mu_lab[in] The desired mu in the lab frame for the emitted particle.
a_energy_out[in] The energy of the emitted outgoing particle.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_cumulative_weight[in] The cumulative multiplicity.
a_checkOrphanProducts[in] If true, associated orphan products are also sampled.
Returns
The weight that the particle is emitted into mu a_mu_lab.

Definition at line 3326 of file MCGIDI_headerSource.hpp.

3327 {
3328
3329 double cumulative_weight1 = 0.0;
3330 if( a_cumulative_weight == nullptr ) a_cumulative_weight = &cumulative_weight1;
3331 double weight1 = 0.0;
3332
3333#ifdef MCGIDI_USE_OUTPUT_CHANNEL
3334 m_outputChannel->angleBiasingViaIntid( this, a_intid, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3335#else
3336
3337 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
3338 (*productIter)->angleBiasingViaIntid( this, a_intid, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3339 }
3340
3341 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_intid == PoPI::Intids::neutron ) ) {
3342 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
3343 DelayedNeutron const *delayedNeutron1 = m_delayedNeutrons[i1];
3344 Product const &product = delayedNeutron1->product( );
3345
3346 product.angleBiasingViaIntid( this, a_intid, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng, *a_cumulative_weight );
3347 }
3348 }
3349#endif
3350
3351 if( a_checkOrphanProducts ) {
3352 for( auto productIter = m_associatedOrphanProducts.begin( ); productIter != m_associatedOrphanProducts.end( ); ++productIter ) {
3353 (*productIter)->angleBiasingViaIntid( this, a_intid, a_temperature, a_energy_in, a_mu_lab, weight1, a_energy_out, a_rng,
3354 *a_cumulative_weight );
3355 }
3356 }
3357
3358 return( weight1 );
3359}
static int constexpr neutron
Definition PoPI.hpp:177

◆ associatedOrphanProductIndices()

LUPI_HOST_DEVICE Vector< std::size_t > MCGIDI::Reaction::associatedOrphanProductIndices ( ) const
inline

Returns the value of the m_associatedOrphanProductIndicex member.

Definition at line 1439 of file MCGIDI.hpp.

◆ crossSection() [1/2]

LUPI_HOST_DEVICE double MCGIDI::Reaction::crossSection ( URR_protareInfos const & a_URR_protareInfos,
double a_temperature,
double a_energy_in ) const

Returns the reaction's cross section for target temperature a_temperature and projectile energy a_energy_in.

Parameters
a_URR_protareInfos[in] URR information.
a_temperature[in] The temperature of the target.
a_energy_in[in] The energy of the projectile.

Definition at line 232 of file MCGIDI_reaction.cc.

232 {
233
234 return( m_protareSingle->reactionCrossSection( m_reactionIndex, a_URR_protareInfos, a_temperature, a_energy_in ) );
235}

◆ crossSection() [2/2]

LUPI_HOST_DEVICE double MCGIDI::Reaction::crossSection ( URR_protareInfos const & a_URR_protareInfos,
std::size_t a_hashIndex,
double a_temperature,
double a_energy_in ) const

Returns the reaction's cross section for target temperature a_temperature and projectile energy a_energy_in.

Parameters
a_URR_protareInfos[in] URR information.
a_hashIndex[in] Specifies the continuous energy or multi-group index.
a_temperature[in] The temperature of the target.
a_energy_in[in] The energy of the projectile.

Definition at line 219 of file MCGIDI_reaction.cc.

219 {
220
221 return( m_protareSingle->reactionCrossSection( m_reactionIndex, a_URR_protareInfos, a_hashIndex, a_temperature, a_energy_in, false ) );
222}

◆ crossSectionAsGIDI_XYs1d()

LUPI_HOST GIDI::Functions::XYs1d MCGIDI::Reaction::crossSectionAsGIDI_XYs1d ( double a_temperature) const

Returns the reaction's cross section as a pointer to a GIDI::Functions::XYs1d instance.

Returns
A GIDI::Functions::XYs1d instance.

Definition at line 243 of file MCGIDI_reaction.cc.

243 {
244
245 return( m_protareSingle->heatedCrossSections( ).reactionCrossSectionAsGIDI_XYs1d( m_reactionIndex, a_temperature ) );
246}

◆ crossSectionThreshold()

LUPI_HOST_DEVICE double MCGIDI::Reaction::crossSectionThreshold ( ) const
inline

Returns the value of the m_crossSectionThreshold.

Definition at line 1411 of file MCGIDI.hpp.

◆ ENDF_MT()

LUPI_HOST_DEVICE int MCGIDI::Reaction::ENDF_MT ( ) const
inline

Returns the value of the m_ENDF_MT.

Definition at line 1403 of file MCGIDI.hpp.

◆ ENDL_C()

LUPI_HOST_DEVICE int MCGIDI::Reaction::ENDL_C ( ) const
inline

Returns the value of the m_ENDL_C.

Definition at line 1404 of file MCGIDI.hpp.

◆ ENDL_S()

LUPI_HOST_DEVICE int MCGIDI::Reaction::ENDL_S ( ) const
inline

Returns the value of the m_ENDL_S.

Definition at line 1405 of file MCGIDI.hpp.

◆ finalQ()

LUPI_HOST_DEVICE double MCGIDI::Reaction::finalQ ( double a_energy) const

Returns the Q-value for projectile energy a_energy.

Parameters
a_URR_protareInfos[in] URR information.
a_hashIndex[in] Specifies the continuous energy or multi-group index.
a_temperature[in] The temperature of the target.
a_energy_in[in] The energy of the projectile.

Definition at line 198 of file MCGIDI_reaction.cc.

198 {
199
200#ifdef MCGIDI_USE_OUTPUT_CHANNEL
201 return( m_outputChannel->finalQ( a_energy ) );
202#else
203 double Q = 0.0;
204 for( auto Q_iter = m_Qs.begin( ); Q_iter != m_Qs.end( ); ++Q_iter ) Q += (*Q_iter)->evaluate( a_energy );
205
206 return( Q );
207#endif
208}

◆ GIDI_reactionIndex()

LUPI_HOST_DEVICE std::size_t MCGIDI::Reaction::GIDI_reactionIndex ( ) const
inline

Returns the value of the m_GIDI_reactionIndex member.

Definition at line 1401 of file MCGIDI.hpp.

◆ hasFission()

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

Returns the value of the m_hasFission.

Definition at line 1408 of file MCGIDI.hpp.

◆ initialStateIndex()

LUPI_HOST_DEVICE int MCGIDI::Reaction::initialStateIndex ( ) const
inline

Returns the value of the m_initialStateIndex member.

Definition at line 1406 of file MCGIDI.hpp.

◆ label()

LUPI_HOST_DEVICE String const & MCGIDI::Reaction::label ( ) const
inline

Returns the value of the m_label.

Definition at line 1402 of file MCGIDI.hpp.

◆ numberOfProducts()

LUPI_HOST_DEVICE std::size_t MCGIDI::Reaction::numberOfProducts ( ) const
inline

Returns the number of products in the m_products member.

Definition at line 1421 of file MCGIDI.hpp.

◆ product()

LUPI_HOST_DEVICE Product const * MCGIDI::Reaction::product ( std::size_t a_index) const
inline

Definition at line 1422 of file MCGIDI.hpp.

1422{ return( m_products[a_index] ); }

Referenced by angleBiasing(), angleBiasingViaIntid(), and sampleProducts().

◆ productAverageMultiplicity()

LUPI_HOST_DEVICE double MCGIDI::Reaction::productAverageMultiplicity ( int a_index,
double a_projectileEnergy ) const

Returns the energy dependent multiplicity for outgoing particle with pops index a_index. The returned value may not be an integer. Energy dependent multiplicity mainly occurs for photons and fission neutrons.

Parameters
a_index[in] The PoPs index of the requested particle.
a_projectileEnergy[in] The energy of the projectile.
Returns
The multiplicity value for the requested particle.

Definition at line 297 of file MCGIDI_reaction.cc.

297 {
298
299 double multiplicity = 0.0;
300
301 if( m_crossSectionThreshold > a_projectileEnergy ) return( multiplicity );
302
303 std::size_t i1 = 0;
304 for( Vector<int>::iterator iter = m_productIndices.begin( ); iter != m_productIndices.end( ); ++iter, ++i1 ) {
305 if( *iter == a_index ) {
306 multiplicity = m_productMultiplicities[i1];
307 break;
308 }
309 }
310
311 if( multiplicity < 0 ) {
312#ifdef MCGIDI_USE_OUTPUT_CHANNEL
313 multiplicity = m_outputChannel->productAverageMultiplicity( a_index, a_projectileEnergy );
314#else
315 multiplicity = 0.0;
316 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
317 multiplicity += (*productIter)->productAverageMultiplicity( a_index, a_projectileEnergy );
318 }
319
320 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_index == m_neutronIndex ) ) {
321 multiplicity += m_totalDelayedNeutronMultiplicity->evaluate( a_projectileEnergy );
322 }
323#endif
324 }
325
326 return( multiplicity );
327}

◆ productAverageMultiplicityViaIntid()

LUPI_HOST_DEVICE double MCGIDI::Reaction::productAverageMultiplicityViaIntid ( int a_intid,
double a_projectileEnergy ) const

Returns the energy dependent multiplicity for outgoing particle with pops intid a_intid. The returned value may not be an integer. Energy dependent multiplicity mainly occurs for photons and fission neutrons.

Parameters
a_intid[in] The PoPs intid of the requested particle.
a_projectileEnergy[in] The energy of the projectile.
Returns
The multiplicity value for the requested particle.

Definition at line 339 of file MCGIDI_reaction.cc.

339 {
340
341 double multiplicity = 0.0;
342
343 if( m_crossSectionThreshold > a_projectileEnergy ) return( multiplicity );
344
345 std::size_t i1 = 0;
346 for( Vector<int>::iterator iter = m_productIntids.begin( ); iter != m_productIntids.end( ); ++iter, ++i1 ) {
347 if( *iter == a_intid ) {
348 multiplicity = m_productMultiplicities[i1];
349 break;
350 }
351 }
352
353 if( multiplicity < 0 ) {
354#ifdef MCGIDI_USE_OUTPUT_CHANNEL
355 multiplicity = m_outputChannel->productAverageMultiplicityViaIntid( a_intid, a_projectileEnergy );
356#else
357 multiplicity = 0.0;
358 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
359 multiplicity += (*productIter)->productAverageMultiplicityViaIntid( a_intid, a_projectileEnergy );
360 }
361
362 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_intid == PoPI::Intids::neutron ) ) {
363 multiplicity += m_totalDelayedNeutronMultiplicity->evaluate( a_projectileEnergy );
364 }
365#endif
366 }
367
368 return( multiplicity );
369}

◆ productIndices()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::productIndices ( ) const
inline

Returns a const reference to the m_productIntids member.

Definition at line 1419 of file MCGIDI.hpp.

◆ productIndicesTransportable()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::productIndicesTransportable ( ) const
inline

Returns a const reference to the m_productIndicesTransportable member.

Definition at line 1432 of file MCGIDI.hpp.

◆ productIntids()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::productIntids ( ) const
inline

Definition at line 1418 of file MCGIDI.hpp.

1418{ return( m_productIntids ); }

◆ productIntidsTransportable()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::productIntidsTransportable ( ) const
inline

Returns a const reference to the m_productIntidsTransportable member.

Definition at line 1430 of file MCGIDI.hpp.

◆ productMultiplicities()

LUPI_HOST_DEVICE int MCGIDI::Reaction::productMultiplicities ( int a_index) const
inline

This method is deprecated. Please use productMultiplicity instead.

Definition at line 1425 of file MCGIDI.hpp.

◆ productMultiplicity()

LUPI_HOST_DEVICE int MCGIDI::Reaction::productMultiplicity ( int a_index) const

Returns the multiplicity for outgoing particle with pops index a_index. If the multiplicity is energy dependent, the returned value is -1. For energy dependent multiplicities it is better to use the method productAverageMultiplicity.

Parameters
a_index[in] The PoPs index of the requested particle.
Returns
The multiplicity value for the requested particle.

Definition at line 257 of file MCGIDI_reaction.cc.

257 {
258
259 std::size_t i1 = 0;
260
261 for( Vector<int>::iterator iter = m_productIndices.begin( ); iter != m_productIndices.end( ); ++iter, ++i1 ) {
262 if( *iter == a_index ) return( m_productMultiplicities[i1] );
263 }
264
265 return( 0 );
266}

Referenced by productMultiplicities().

◆ productMultiplicityViaIntid()

LUPI_HOST_DEVICE int MCGIDI::Reaction::productMultiplicityViaIntid ( int a_intid) const

Returns the multiplicity for outgoing particle with pops intid a_intid. If the multiplicity is energy dependent, the returned value is -1. For energy dependent multiplicities it is better to use the method productAverageMultiplicity.

Parameters
a_intid[in] The PoPs intid of the requested particle.
Returns
The multiplicity value for the requested particle.

Definition at line 276 of file MCGIDI_reaction.cc.

276 {
277
278 std::size_t i1 = 0;
279
280 for( Vector<int>::iterator iter = m_productIntids.begin( ); iter != m_productIntids.end( ); ++iter, ++i1 ) {
281 if( *iter == a_intid ) return( m_productMultiplicities[i1] );
282 }
283
284 return( 0 );
285}

◆ projectileMass()

LUPI_HOST_DEVICE double MCGIDI::Reaction::projectileMass ( ) const
inline

Returns the value of the m_projectileMass.

Definition at line 1409 of file MCGIDI.hpp.

◆ protareSingle()

◆ reactionIndex()

◆ sampleNullProducts()

template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void MCGIDI::Reaction::sampleNullProducts ( Protare const & a_protare,
double a_projectileEnergy,
Sampling::Input & a_input,
RNG && a_rng,
PUSHBACK && a_push_back,
Sampling::ProductHandler & a_products )
inlinestatic

This method adds a null product to a_products. When running in multi-group mode, a sampled reaction may be rejected if the threshold is in the multi-group that the projectile is in. If this happens, only null products should be returned. This type of behavior was need in MCAPM but is probably not needed for MCGIDI.

Parameters
a_protare[in] The Protare this Reaction belongs to.
a_projectileEnergy[in] The energy of the projectile.
a_input[in] Sample options requested by user.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_products[in] The object to add all sampled products to.

Definition at line 3240 of file MCGIDI_headerSource.hpp.

3241 {
3242
3243 a_input.m_GRIN_intermediateResidual = -1;
3244 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
3245 a_input.m_dataInTargetFrame = false;
3246 a_input.m_frame = GIDI::Frame::lab;
3247 a_input.m_delayedNeutronIndex = -1;
3248 a_input.m_delayedNeutronDecayRate = 0.0;
3249
3250 a_input.m_energyOut1 = a_projectileEnergy;
3251 a_input.m_mu = 1.0;
3252 a_input.m_phi = 0.0;
3253
3254 a_products.add( a_projectileEnergy, a_protare.projectileIntid( ), a_protare.projectileIndex( ), a_protare.projectileUserIndex( ), a_protare.projectileMass( ), a_input, a_rng, a_push_back, false );
3255}

◆ sampleProducts()

template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void MCGIDI::Reaction::sampleProducts ( Protare const * a_protare,
Sampling::Input & a_input,
RNG && a_rng,
PUSHBACK && a_push_back,
Sampling::ProductHandler & a_products,
bool a_checkOrphanProducts = true ) const
inline

This method adds sampled products to a_products.

Parameters
a_protare[in] The Protare this Reaction belongs to.
a_input[in] Sample options requested by user.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_products[in] The object to add all sampled products to.
a_checkOrphanProducts[in] If true, associated orphan products are also sampled.

Definition at line 2966 of file MCGIDI_headerSource.hpp.

2967 {
2968
2969 double projectileEnergy = a_input.modelEnergy( );
2970
2971 a_input.m_GRIN_intermediateResidual = -1;
2972 a_input.m_reaction = this;
2973 a_input.m_projectileMass = m_projectileMass;
2974 a_input.m_targetMass = m_targetMass;
2975 a_input.m_relativeBeta = MCGIDI_particleBeta( m_projectileMass, projectileEnergy );
2976
2977 if( m_GRIN_specialSampleProducts ) {
2978 if( m_GRIN_capture != nullptr ) {
2979 if( projectileEnergy < m_GRIN_maximumCaptureIncidentEnergy ) {
2980 if( m_GRIN_capture->sampleProducts( (ProtareSingle const *) a_protare, projectileEnergy, a_input, a_rng, a_push_back, a_products ) ) {
2981 return;
2982 }
2983 } }
2984 else if( m_GRIN_inelastic != nullptr ) {
2985 if( m_GRIN_inelastic->sampleProducts( (ProtareSingle const *) a_protare, projectileEnergy, a_input, a_rng, a_push_back, a_products ) ) {
2986 return;
2987 }
2988 }
2989 }
2990
2991#ifdef MCGIDI_USE_OUTPUT_CHANNEL
2992 m_outputChannel->sampleProducts( m_protareSingle, projectileEnergy, a_input, a_rng, a_push_back, a_products );
2993#else
2994 if( m_hasFinalStatePhotons ) {
2995 double random = a_rng( );
2996 double cumulative = 0.0;
2997 bool sampled = false;
2998 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
2999 cumulative += (*productIter)->multiplicity( )->evaluate( projectileEnergy );
3000 if( cumulative >= random ) {
3001 (*productIter)->sampleFinalState( m_protareSingle, projectileEnergy, a_input, a_rng, a_push_back, a_products );
3002 sampled = true;
3003 break;
3004 }
3005 }
3006 if( !sampled ) { // BRB: FIXME: still need to code for continuum photon.
3007 } }
3008 else {
3009 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
3010 (*productIter)->sampleProducts( m_protareSingle, projectileEnergy, a_input, a_rng, a_push_back, a_products );
3011 }
3012 }
3013
3014 if( m_totalDelayedNeutronMultiplicity != nullptr ) {
3015 double totalDelayedNeutronMultiplicity = m_totalDelayedNeutronMultiplicity->evaluate( projectileEnergy );
3016
3017 if( a_rng( ) < totalDelayedNeutronMultiplicity ) { // Assumes that totalDelayedNeutronMultiplicity < 1.0, which it is.
3018 double sum = 0.0;
3019
3020 totalDelayedNeutronMultiplicity *= a_rng( );
3021 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
3022 DelayedNeutron const *delayedNeutron1 = m_delayedNeutrons[i1];
3023 Product const &product = delayedNeutron1->product( );
3024
3025 sum += product.multiplicity( )->evaluate( projectileEnergy );
3026 if( sum >= totalDelayedNeutronMultiplicity ) {
3027 product.distribution( )->sample( projectileEnergy, a_input, a_rng );
3028 a_input.m_delayedNeutronIndex = delayedNeutron1->delayedNeutronIndex( );
3029 a_input.m_delayedNeutronDecayRate = delayedNeutron1->rate( );
3030 a_products.add( a_input.energy( ), product.intid( ), product.index( ), product.userParticleIndex( ), product.mass( ), a_input, a_rng, a_push_back, false );
3031 break;
3032 }
3033 }
3034 }
3035 }
3036
3037 if( m_fissionResiduaIntid != -1 ) { // Special treatment to add 2 ENDL 99120 or 99125 products.
3038 a_input.setSampledType( MCGIDI::Sampling::SampledType::unspecified );
3039 a_input.m_frame = GIDI::Frame::lab;
3040 a_input.m_energyOut1 = 0.0;
3041 a_input.m_mu = 0.0;
3042 a_input.m_phi = 0.0;
3043 a_input.m_delayedNeutronIndex = -1;
3044 a_input.m_delayedNeutronDecayRate = 0.0;
3045 a_products.add( 0.0, m_fissionResiduaIntid, m_fissionResiduaIndex, m_fissionResiduaUserIndex, m_fissionResidualMass, a_input, a_rng, a_push_back, false );
3046 a_products.add( 0.0, m_fissionResiduaIntid, m_fissionResiduaIndex, m_fissionResiduaUserIndex, m_fissionResidualMass, a_input, a_rng, a_push_back, false );
3047 }
3048
3049#endif
3050
3051 if( a_checkOrphanProducts ) {
3052 for( auto productIter = m_associatedOrphanProducts.begin( ); productIter != m_associatedOrphanProducts.end( ); ++productIter ) {
3053 (*productIter)->sampleProducts( m_protareSingle, projectileEnergy, a_input, a_rng, a_push_back, a_products );
3054 }
3055 }
3056}
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71

◆ serialize()

LUPI_HOST_DEVICE void MCGIDI::Reaction::serialize ( 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 530 of file MCGIDI_reaction.cc.

530 {
531
532 DATA_MEMBER_SIZE_T( m_GIDI_reactionIndex, a_buffer, a_mode );
533 DATA_MEMBER_STRING( m_label, a_buffer, a_mode );
534 DATA_MEMBER_INT( m_ENDF_MT, a_buffer, a_mode );
535 DATA_MEMBER_INT( m_ENDL_C, a_buffer, a_mode );
536 DATA_MEMBER_INT( m_ENDL_S, a_buffer, a_mode );
537 DATA_MEMBER_INT( m_initialStateIndex, a_buffer, a_mode );
538 DATA_MEMBER_INT( m_neutronIndex, a_buffer, a_mode );
539 DATA_MEMBER_CAST( m_hasFission, a_buffer, a_mode, bool );
540 DATA_MEMBER_DOUBLE( m_projectileMass, a_buffer, a_mode );
541 DATA_MEMBER_DOUBLE( m_targetMass, a_buffer, a_mode );
542 DATA_MEMBER_DOUBLE( m_crossSectionThreshold, a_buffer, a_mode );
543 DATA_MEMBER_DOUBLE( m_twoBodyThreshold, a_buffer, a_mode );
544 DATA_MEMBER_CAST( m_hasFinalStatePhotons, a_buffer, a_mode, bool );
545 DATA_MEMBER_INT( m_fissionResiduaIntid, a_buffer, a_mode );
546 DATA_MEMBER_INT( m_fissionResiduaIndex, a_buffer, a_mode );
547 DATA_MEMBER_INT( m_fissionResiduaUserIndex, a_buffer, a_mode );
548 serializeFissionResiduals( m_fissionResiduals, a_buffer, a_mode );
549 DATA_MEMBER_DOUBLE( m_fissionResidualMass, a_buffer, a_mode );
550
551 DATA_MEMBER_VECTOR_INT( m_productIntids, a_buffer, a_mode );
552 DATA_MEMBER_VECTOR_INT( m_productIndices, a_buffer, a_mode );
553 DATA_MEMBER_VECTOR_INT( m_userProductIndices, a_buffer, a_mode );
554 DATA_MEMBER_VECTOR_INT( m_productMultiplicities, a_buffer, a_mode );
555 DATA_MEMBER_VECTOR_INT( m_productIntidsTransportable, a_buffer, a_mode );
556 DATA_MEMBER_VECTOR_INT( m_productIndicesTransportable, a_buffer, a_mode );
557 DATA_MEMBER_VECTOR_INT( m_userProductIndicesTransportable, a_buffer, a_mode );
558
559#ifdef MCGIDI_USE_OUTPUT_CHANNEL
560 bool haveChannel = m_outputChannel != nullptr;
561 DATA_MEMBER_CAST( haveChannel, a_buffer, a_mode, bool );
562 if( haveChannel ) {
563 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
564 if (a_buffer.m_placement != nullptr) {
565 m_outputChannel = new(a_buffer.m_placement) OutputChannel();
566 a_buffer.incrementPlacement( sizeof(OutputChannel));
567 }
568 else {
569 m_outputChannel = new OutputChannel();
570 } }
571 else if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
572 a_buffer.incrementPlacement( sizeof( OutputChannel ) );
573 }
574 m_outputChannel->serialize( a_buffer, a_mode );
575 }
576#else
577 serializeQs( a_buffer, a_mode, m_Qs );
578 serializeProducts( a_buffer, a_mode, m_products );
579 m_totalDelayedNeutronMultiplicity = serializeFunction1d( a_buffer, a_mode, m_totalDelayedNeutronMultiplicity );
580 serializeDelayedNeutrons( a_buffer, a_mode, m_delayedNeutrons );
581#endif
582
583 DATA_MEMBER_VECTOR_SIZE_T( m_associatedOrphanProductIndices, a_buffer, a_mode );
584
585 std::size_t vectorSize = m_associatedOrphanProducts.size( );
586 int vectorSizeInt = (int) vectorSize;
587 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
588 vectorSize = (std::size_t) vectorSizeInt;
589 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
590 m_associatedOrphanProducts.reserve( vectorSize, &a_buffer.m_placement ); }
591 else if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
592 a_buffer.m_placement += m_associatedOrphanProducts.internalSize( );
593 }
594
595 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
596 m_protareSingle = nullptr;
597 m_reactionIndex = MCGIDI_nullReaction;
598 }
599
600 DATA_MEMBER_CAST( m_GRIN_specialSampleProducts, a_buffer, a_mode, bool );
601 DATA_MEMBER_DOUBLE( m_GRIN_inelasticThreshold, a_buffer, a_mode );
602 DATA_MEMBER_DOUBLE( m_GRIN_maximumCaptureIncidentEnergy, a_buffer, a_mode );
603
604 if( m_GRIN_specialSampleProducts ) {
605 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
606 if( a_buffer.m_placement != nullptr ) {
607 if( m_ENDF_MT == 91 ) {
608 m_GRIN_inelastic = new(a_buffer.m_placement) GRIN_inelastic;
609 a_buffer.incrementPlacement( sizeof( GRIN_inelastic ) ); }
610 else {
611 m_GRIN_capture = new(a_buffer.m_placement) GRIN_capture;
612 a_buffer.incrementPlacement( sizeof( GRIN_capture ) );
613 } }
614 else {
615 if( m_ENDF_MT == 91 ) {
616 m_GRIN_inelastic = new GRIN_inelastic( ); }
617 else {
618 m_GRIN_capture = new GRIN_capture( );
619 }
620 }
621 }
622
623 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
624 if( m_ENDF_MT == 91 ) {
625 a_buffer.incrementPlacement( sizeof( GRIN_inelastic ) ); }
626 else {
627 a_buffer.incrementPlacement( sizeof( GRIN_capture ) );
628 }
629 }
630
631 if( m_ENDF_MT == 91 ) {
632 m_GRIN_inelastic->serialize( a_buffer, a_mode ); }
633 else {
634 m_GRIN_capture->serialize( a_buffer, a_mode );
635 }
636 }
637}
#define DATA_MEMBER_STRING(member, buf, mode)
#define DATA_MEMBER_VECTOR_INT(member, buf, mode)
#define DATA_MEMBER_CAST(member, buf, mode, someType)
#define DATA_MEMBER_SIZE_T(member, buf, mode)
#define DATA_MEMBER_VECTOR_SIZE_T(member, buf, mode)
#define DATA_MEMBER_DOUBLE(member, buf, mode)
#define DATA_MEMBER_INT( member, buf, mode)
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
LUPI_HOST_DEVICE void serializeFissionResiduals(GIDI::Construction::FissionResiduals &a_fissionResiduals, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serializeDelayedNeutrons(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< DelayedNeutron * > &a_delayedNeutrons)
LUPI_HOST_DEVICE void serializeQs(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< Functions::Function1d_d1 * > &a_Qs)
LUPI_HOST_DEVICE Functions::Function1d * serializeFunction1d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d *a_function1d)
LUPI_HOST_DEVICE void serializeProducts(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< Product * > &a_products)

◆ setModelDBRC_data()

LUPI_HOST void MCGIDI::Reaction::setModelDBRC_data ( Sampling::Upscatter::ModelDBRC_data * a_modelDBRC_data)

This method calls the **setModelDBRC_data* method on the first product of this with a_modelDBRC_data.

Parameters
a_modelDBRC_data[in] The instance storing data needed to treat the DRRC upscatter mode.

Definition at line 445 of file MCGIDI_reaction.cc.

445 {
446
447#ifdef MCGIDI_USE_OUTPUT_CHANNEL
448 m_outputChannel->setModelDBRC_data( a_modelDBRC_data );
449#else
450 m_products[0]->setModelDBRC_data( a_modelDBRC_data );
451#endif
452}

◆ setOrphanProductData()

LUPI_HOST void MCGIDI::Reaction::setOrphanProductData ( std::vector< std::size_t > const & a_associatedOrphanProductIndcies,
std::vector< Product * > const & a_associatedOrphanProducts )

Adds the contents of a_associatedOrphanProductIndcies and a_associatedOrphanProducts to this instance.

Parameters
a_associatedOrphanProductIndcies[in] The list of indices of the orphanProduct reaction that make up the product in a_associatedOrphanProducts.
a_associatedOrphanProducts[in] The list of pointers to the associated orphan products.

Definition at line 510 of file MCGIDI_reaction.cc.

511 {
512
513 m_associatedOrphanProductIndices.reserve( a_associatedOrphanProductIndcies.size( ) );
514 for( auto iter = a_associatedOrphanProductIndcies.begin( ); iter != a_associatedOrphanProductIndcies.end( ); ++iter )
515 m_associatedOrphanProductIndices.push_back( *iter );
516
517 m_associatedOrphanProducts.reserve( a_associatedOrphanProducts.size( ) );
518 for( auto iter = a_associatedOrphanProducts.begin( ); iter != a_associatedOrphanProducts.end( ); ++iter )
519 m_associatedOrphanProducts.push_back( *iter );
520}

◆ setUserParticleIndex()

LUPI_HOST void MCGIDI::Reaction::setUserParticleIndex ( 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 378 of file MCGIDI_reaction.cc.

378 {
379
380#ifdef MCGIDI_USE_OUTPUT_CHANNEL
381 m_outputChannel->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
382#else
383 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
384 (*productIter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
385 }
386
387 for( auto iter = m_delayedNeutrons.begin( ); iter != m_delayedNeutrons.end( ); ++iter )
388 (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
389#endif
390
391 for( std::size_t i1 = 0; i1 < m_productIndices.size( ); ++i1 ) {
392 if( m_productIndices[i1] == a_particleIndex ) m_userProductIndices[i1] = a_userParticleIndex;
393 }
394
395 for( std::size_t i1 = 0; i1 < m_productIndicesTransportable.size( ); ++i1 ) {
396 if( m_productIndicesTransportable[i1] == a_particleIndex ) m_userProductIndicesTransportable[i1] = a_userParticleIndex;
397 }
398
399 if( a_particleIndex == m_fissionResiduaIndex ) m_fissionResiduaUserIndex = a_userParticleIndex;
400
401 if( m_GRIN_capture != nullptr ) m_GRIN_capture->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
402 if( m_GRIN_inelastic != nullptr ) m_GRIN_inelastic->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
403}

◆ setUserParticleIndexViaIntid()

LUPI_HOST void MCGIDI::Reaction::setUserParticleIndexViaIntid ( int a_particleIntid,
int a_userParticleIndex )

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

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

Definition at line 412 of file MCGIDI_reaction.cc.

412 {
413
414#ifdef MCGIDI_USE_OUTPUT_CHANNEL
415 m_outputChannel->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
416#else
417 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
418 (*productIter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
419 }
420
421 for( auto iter = m_delayedNeutrons.begin( ); iter != m_delayedNeutrons.end( ); ++iter )
422 (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
423#endif
424
425 for( std::size_t i1 = 0; i1 < m_productIntids.size( ); ++i1 ) {
426 if( m_productIntids[i1] == a_particleIntid ) m_userProductIndices[i1] = a_userParticleIndex;
427 }
428
429 for( std::size_t i1 = 0; i1 < m_productIntidsTransportable.size( ); ++i1 ) {
430 if( m_productIntidsTransportable[i1] == a_particleIntid ) m_userProductIndicesTransportable[i1] = a_userParticleIndex;
431 }
432
433 if( a_particleIntid == m_fissionResiduaIntid ) m_fissionResiduaUserIndex = a_userParticleIndex;
434
435 if( m_GRIN_capture != nullptr ) m_GRIN_capture->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
436 if( m_GRIN_inelastic != nullptr ) m_GRIN_inelastic->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
437}

◆ targetMass()

LUPI_HOST_DEVICE double MCGIDI::Reaction::targetMass ( ) const
inline

Returns the value of the m_targetMass.

Definition at line 1410 of file MCGIDI.hpp.

◆ twoBodyThreshold()

LUPI_HOST_DEVICE double MCGIDI::Reaction::twoBodyThreshold ( ) const
inline

Returns the value of the m_twoBodyThreshold member.

Definition at line 1412 of file MCGIDI.hpp.

◆ updateProtareSingleInfo()

LUPI_HOST_DEVICE void MCGIDI::Reaction::updateProtareSingleInfo ( ProtareSingle * a_protareSingle,
std::size_t a_reactionIndex )
inline

Definition at line 1395 of file MCGIDI.hpp.

1395 {
1396 m_protareSingle = a_protareSingle;
1397 m_reactionIndex = a_reactionIndex;
1398 }

◆ userProductIndices()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::userProductIndices ( ) const
inline

Returns a const reference to the m_productIndices member.

Definition at line 1420 of file MCGIDI.hpp.

◆ userProductIndicesTransportable()

LUPI_HOST_DEVICE Vector< int > const & MCGIDI::Reaction::userProductIndicesTransportable ( ) const
inline

Definition at line 1434 of file MCGIDI.hpp.

1434{ return( m_userProductIndicesTransportable ); }

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