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

#include <MCGIDI.hpp>

Inheritance diagram for MCGIDI::ProtareSingle:

Public Member Functions

LUPI_HOST_DEVICE ProtareSingle ()
LUPI_HOST ProtareSingle (LUPI::StatusMessageReporting &a_smr, GIDI::ProtareSingle 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 ~ProtareSingle ()
LUPI_HOST_DEVICE bool isPhotoAtomic () const
LUPI_HOST_DEVICE bool continuousEnergy () const
LUPI_HOST_DEVICE bool fixedGrid () const
LUPI_HOST_DEVICE HeatedCrossSectionsContinuousEnergy const & heatedCrossSections () const
LUPI_HOST_DEVICE HeatedCrossSectionsContinuousEnergyheatedCrossSections ()
LUPI_HOST_DEVICE HeatedCrossSectionsMultiGroup const & heatedMultigroupCrossSections () const
LUPI_HOST_DEVICE HeatedCrossSectionsMultiGroupheatedMultigroupCrossSections ()
LUPI_HOST_DEVICE const Vector< NuclideGammaBranchStateInfo * > & nuclideGammaBranchStateInfos () const
LUPI_HOST_DEVICE const Vector< NuclideGammaBranchInfo * > & branches () const
LUPI_HOST_DEVICE Vector< Reaction * > const & reactions () const
LUPI_HOST_DEVICE Vector< Reaction * > const & orphanProducts () const
template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void sampleBranchingGammas (Sampling::Input &a_input, double a_projectileEnergy, int a_initialStateIndex, RNG &&a_rng, PUSHBACK &&push_back, Sampling::ProductHandler &a_products) 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 String interaction () const
LUPI_HOST_DEVICE bool hasIncoherentDoppler () const
LUPI_HOST_DEVICE int URR_index () const
LUPI_HOST_DEVICE void setURR_index (int a_URR_index)
LUPI_HOST_DEVICE bool inURR (double a_energy) 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 bool sampleTargetBetaForUpscatterModelA (Sampling::Input &a_input, RNG &&a_rng) 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 bool upscatterModelASupported () const
LUPI_HOST_DEVICE Vector< double > const & upscatterModelAGroupEnergies () const
LUPI_HOST_DEVICE Vector< double > const & upscatterModelAGroupVelocities () const
LUPI_HOST_DEVICE Vector< double > const & upscatterModelACrossSection () 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 1609 of file MCGIDI.hpp.

Constructor & Destructor Documentation

◆ ProtareSingle() [1/2]

LUPI_HOST_DEVICE MCGIDI::ProtareSingle::ProtareSingle ( )

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

Definition at line 1272 of file MCGIDI_protare.cc.

1272 :
1274 m_URR_index( -1 ),
1275 m_hasURR_probabilityTables( false ),
1276 m_URR_domainMin( -1.0 ),
1277 m_URR_domainMax( -1.0 ),
1278 m_upscatterModelASupported( false ),
1279 m_projectileMultiGroupBoundaries( 0 ),
1280 m_projectileMultiGroupBoundariesCollapsed( 0 ),
1281 m_reactions( 0 ),
1282 m_orphanProducts( 0 ) {
1283
1284}
LUPI_HOST_DEVICE Protare(ProtareType a_protareType)

Referenced by protare(), protare(), and protareWithReaction().

◆ ProtareSingle() [2/2]

LUPI_HOST MCGIDI::ProtareSingle::ProtareSingle ( LUPI::StatusMessageReporting & a_smr,
GIDI::ProtareSingle 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 )
Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_protare[in] The GIDI::Protare whose data is to be used to construct this.
a_pops[in] A PoPs Database instance used to get particle intids and possibly other particle information.
a_settings[in] Used to pass user options to the this to instruct it which data are desired.
a_particles[in] List of transporting particles and their information (e.g., multi-group boundaries and fluxes).
a_domainHash[in] The hash data used when looking up a cross section.
a_temperatureInfos[in] The list of temperature data to extract from a_protare.
a_reactionsToExclude[in] A list of reaction to not include in the MCGIDI::Protare.
a_reactionsToExcludeOffset[in] The starting index for the reactions in this ProtareSingle.
a_allowFixedGrid[in] For internal (i.e., MCGIDI) use only. Users must use the default value.

Definition at line 1299 of file MCGIDI_protare.cc.

1302 :
1303 Protare( ProtareType::single, a_protare, a_settings, a_pops ),
1304 m_interaction( a_protare.interaction( ).c_str( ) ),
1305 m_URR_index( -1 ),
1306 m_hasURR_probabilityTables( false ),
1307 m_URR_domainMin( -1.0 ),
1308 m_URR_domainMax( -1.0 ),
1309 m_domainHash( a_domainHash ),
1310 m_upscatterModelASupported( ( projectileIntid( ) != PoPI::Intids::photon ) &&
1313 m_projectileMultiGroupBoundaries( 0 ),
1314 m_projectileMultiGroupBoundariesCollapsed( 0 ),
1315 m_reactions( 0 ),
1316 m_orphanProducts( 0 ),
1317 m_isPhotoAtomic( a_protare.isPhotoAtomic( ) ),
1318 m_heatedCrossSections( ),
1319 m_heatedMultigroupCrossSections( ) {
1320
1321 a_protare.updateReactionIndices( 0 ); // This is not correct as the offset should be passed as an arguent.
1322 PoPI::Database const &pops = a_protare.protare( 0 )->internalPoPs( );
1323
1324 if( !a_protare.isPhotoAtomic( ) ) {
1325 std::set<std::string> incompleteParticles;
1326 a_protare.incompleteParticles( a_settings, incompleteParticles );
1327 for( auto particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( ); ++particle ) {
1328 if( incompleteParticles.count( particle->first ) != 0 ) {
1329 std::string message = "Requested particle '" + particle->first + "' is incomplete in '" + a_protare.realFileName( ) + "'.";
1330 if( a_settings.throwOnError( ) ) {
1331 throw std::runtime_error( message.c_str( ) ); }
1332 else {
1333 smr_setReportError2p( a_smr.smr( ), 0, 0, message.c_str( ) );
1334 }
1335 }
1336 }
1337 }
1338
1339 SetupInfo setupInfo( *this, a_protare, a_pops, pops );
1340 setupInfo.m_formatVersion = a_protare.formatVersion( );
1341 setupInfo.m_GRIN_continuumGammas = a_protare.GRIN_continuumGammas2( );
1342
1343 GIDI::Transporting::Particles particles;
1344 for( std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( ); ++particle ) {
1345 setupInfo.m_particleIntids[particle->first] = MCGIDI_popsIntid( pops, particle->first );
1346 setupInfo.m_particleIndices[particle->first] = MCGIDI_popsIndex( a_pops, particle->first );
1347
1348 if( ( m_interaction == GIDI_MapInteractionAtomicChars ) &&
1349 !( ( particle->first == PoPI::IDs::photon ) || ( particle->first == PoPI::IDs::electron ) ) ) continue;
1350 particles.add( particle->second );
1351 }
1352
1353 GIDI::Transporting::MG multiGroupSettings( a_settings.projectileID( ), GIDI::Transporting::Mode::MonteCarloContinuousEnergy, a_settings.delayedNeutrons( ) );
1354 multiGroupSettings.setThrowOnError( a_settings.throwOnError( ) );
1355
1356 setupInfo.m_distributionLabel = a_temperatureInfos[0].griddedCrossSection( );
1357
1358 a_settings.styles( &a_protare.styles( ) );
1359
1360 switch( a_settings.crossSectionLookupMode( ) ) {
1362 m_continuousEnergy = true;
1363 break;
1365 m_continuousEnergy = false;
1366 if( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::B ) {
1367 multiGroupSettings.setMode( GIDI::Transporting::Mode::multiGroupWithSnElasticUpScatter ); }
1368 else {
1369 multiGroupSettings.setMode( GIDI::Transporting::Mode::multiGroup );
1370 }
1371 break;
1372 default :
1373 throw std::runtime_error( "ProtareSingle::ProtareSingle: invalid lookupMode" );
1374 }
1375 m_fixedGrid = a_allowFixedGrid && ( a_protare.projectile( ).ID( ) == PoPI::IDs::photon ) && ( a_settings.fixedGridPoints( ).size( ) > 0 );
1376
1377 setupNuclideGammaBranchStateInfos( setupInfo, a_protare, a_settings.makePhotonEmissionProbabilitiesOne( ),
1378 a_settings.zeroNuclearLevelEnergyWidth( ) );
1379 convertACE_URR_probabilityTablesFromGIDI( a_protare, a_settings, setupInfo );
1380
1381 if( ( a_settings.crossSectionLookupMode( ) == Transporting::LookupMode::Data1d::multiGroup ) ||
1382 ( a_settings.other1dDataLookupMode( ) == Transporting::LookupMode::Data1d::multiGroup ) ) {
1383
1384 GIDI::Suite const *transportables = nullptr;
1385 if( setupInfo.m_formatVersion.major( ) > 1 ) {
1386 GIDI::Styles::HeatedMultiGroup const &heatedMultiGroup = *a_protare.styles( ).get<GIDI::Styles::HeatedMultiGroup>( a_settings.label( ) );
1387 transportables = &heatedMultiGroup.transportables( ); }
1388 else {
1389 std::vector<GIDI::Suite::const_iterator> tags = a_protare.styles( ).findAllOfMoniker( GIDI_multiGroupStyleChars );
1390
1391 if( tags.size( ) != 1 ) throw std::runtime_error( "MCGIDI::ProtareSingle::ProtareSingle: What is going on here?" );
1392 GIDI::Styles::MultiGroup const &multiGroup = static_cast<GIDI::Styles::MultiGroup const &>( **tags[0] );
1393 transportables = &multiGroup.transportables( );
1394 }
1395
1396 GIDI::Transportable const transportable = *transportables->get<GIDI::Transportable>( a_protare.projectile( ).ID( ) );
1397 m_projectileMultiGroupBoundaries = transportable.groupBoundaries( );
1398 GIDI::Transporting::Particle const *particle = a_particles.particle( a_protare.projectile( ).ID( ) );
1399 m_projectileMultiGroupBoundariesCollapsed = particle->multiGroup( ).boundaries( );
1400 }
1401
1402 std::vector<GIDI::Reaction const *> GIDI_reactions;
1403 std::set<std::string> product_ids;
1404 std::set<std::string> product_ids_transportable;
1405 GIDI::Reaction const *nuclearPlusCoulombInterferenceReaction = nullptr;
1406 if( a_settings.nuclearPlusCoulombInterferenceOnly( ) ) nuclearPlusCoulombInterferenceReaction = a_protare.nuclearPlusCoulombInterferenceOnlyReaction( );
1407
1408 for( std::size_t reactionIndex = 0; reactionIndex < a_protare.reactions( ).size( ); ++reactionIndex ) {
1409 if( a_reactionsToExclude.find( reactionIndex + a_reactionsToExcludeOffset ) != a_reactionsToExclude.end( ) ) continue;
1410
1411 GIDI::Reaction const *GIDI_reaction = a_protare.reaction( reactionIndex );
1412
1413 if( !GIDI_reaction->active( ) ) continue;
1414
1415 if( m_continuousEnergy ) {
1416 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue; }
1417 else {
1418 GIDI::Vector multi_group_cross_section = GIDI_reaction->multiGroupCrossSection( a_smr, multiGroupSettings, a_temperatureInfos[0] );
1419 GIDI::Vector vector = GIDI::collapse( multi_group_cross_section, a_settings, a_particles, 0.0 );
1420
1421 std::size_t i1 = 0;
1422 for( ; i1 < vector.size( ); ++i1 ) if( vector[i1] != 0.0 ) break;
1423 if( i1 == vector.size( ) ) continue;
1424 }
1425 if( a_settings.ignoreENDF_MT5( ) && ( GIDI_reaction->ENDF_MT( ) == 5 ) && ( a_reactionsToExclude.size( ) == 0 ) ) continue;
1426
1427 GIDI_reaction->productIDs( product_ids, particles, false );
1428 GIDI_reaction->productIDs( product_ids_transportable, particles, true );
1429
1430 if( a_settings.nuclearPlusCoulombInterferenceOnly( ) && GIDI_reaction->RutherfordScatteringPresent( ) ) {
1431 if( nuclearPlusCoulombInterferenceReaction != nullptr ) GIDI_reactions.push_back( nuclearPlusCoulombInterferenceReaction ); }
1432 else {
1433 GIDI_reactions.push_back( GIDI_reaction );
1434 }
1435 }
1436
1437 bool zeroReactions = GIDI_reactions.size( ) == 0; // Happens when all reactions are skipped in the prior loop.
1438 if( zeroReactions ) GIDI_reactions.push_back( a_protare.reaction( 0 ) ); // Special case where no reaction in the protare is wanted so the first one is used but its cross section is set to 0.0 at all energies.
1439
1440 setupInfo.m_reactionType = Transporting::Reaction::Type::Reactions;
1441 m_reactions.reserve( GIDI_reactions.size( ) );
1442 for( auto GIDI_reaction = GIDI_reactions.begin( ); GIDI_reaction != GIDI_reactions.end( ); ++GIDI_reaction ) {
1443 setupInfo.m_reaction = *GIDI_reaction;
1444 setupInfo.m_isPairProduction = (*GIDI_reaction)->isPairProduction( );
1445 setupInfo.m_isPhotoAtomicIncoherentScattering = (*GIDI_reaction)->isPhotoAtomicIncoherentScattering( );
1446 setupInfo.m_initialStateIndex = -1;
1447 Reaction *reaction = new Reaction( **GIDI_reaction, setupInfo, a_settings, particles, a_temperatureInfos );
1448 setupInfo.m_initialStateIndices[(*GIDI_reaction)->label( )] = setupInfo.m_initialStateIndex;
1449 reaction->updateProtareSingleInfo( this, m_reactions.size( ) );
1450 m_reactions.push_back( reaction );
1451 }
1452
1453 std::set<int> product_intids;
1454 std::set<int> product_indices;
1455 for( std::set<std::string>::iterator iter = product_ids.begin( ); iter != product_ids.end( ); ++iter ) {
1456 product_intids.insert( MCGIDI_popsIntid( pops, *iter ) );
1457 product_indices.insert( MCGIDI_popsIndex( a_pops, *iter ) );
1458 }
1459 std::set<int> product_intids_transportable;
1460 std::set<int> product_indices_transportable;
1461 for( std::set<std::string>::iterator iter = product_ids_transportable.begin( ); iter != product_ids_transportable.end( ); ++iter ) {
1462 product_intids_transportable.insert( MCGIDI_popsIntid( pops, *iter ) );
1463 product_indices_transportable.insert( MCGIDI_popsIndex( a_pops, *iter ) );
1464 }
1465 productIntidsAndIndices( product_intids, product_intids_transportable, product_indices, product_indices_transportable );
1466
1467 if( a_settings.sampleNonTransportingParticles( ) || particles.hasParticle( PoPI::IDs::photon ) ) {
1468 setupInfo.m_reactionType = Transporting::Reaction::Type::OrphanProducts;
1469 m_orphanProducts.reserve( a_protare.orphanProducts( ).size( ) );
1470 std::vector< std::vector<std::size_t> > associatedOrphanProductIndices( m_reactions.size( ) );
1471
1472 for( std::size_t orphanProductIndex = 0; orphanProductIndex < a_protare.orphanProducts( ).size( ); ++orphanProductIndex ) {
1473 GIDI::Reaction const *GIDI_reaction = a_protare.orphanProduct( orphanProductIndex );
1474
1475 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue;
1476
1477 setupInfo.m_reaction = GIDI_reaction;
1478 Reaction *orphanProductReaction = new Reaction( *GIDI_reaction, setupInfo, a_settings, particles, a_temperatureInfos );
1479 orphanProductReaction->updateProtareSingleInfo( this, m_orphanProducts.size( ) );
1480 m_orphanProducts.push_back( orphanProductReaction );
1481
1482 GIDI::Functions::Reference1d const *reference( GIDI_reaction->crossSection( ).get<GIDI::Functions::Reference1d>( 0 ) );
1483 std::string xlink = reference->xlink( );
1484 GUPI::Ancestry const *ancestry = a_protare.findInAncestry( xlink );
1485 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 1." );
1486 ancestry = ancestry->ancestor( );
1487 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 2." );
1488 if( ancestry->moniker( ) != GIDI_crossSectionSumChars ) {
1489 ancestry = ancestry->ancestor( );
1490 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 3." );
1491 }
1492 GIDI::Sums::CrossSectionSum const *crossSectionSum = static_cast<GIDI::Sums::CrossSectionSum const *>( ancestry );
1493 GIDI::Sums::Summands const &summands = crossSectionSum->summands( );
1494 for( std::size_t i1 = 0; i1 < summands.size( ); ++i1 ) {
1495 GIDI::Sums::Summand::Base const *summand = summands[i1];
1496
1497 ancestry = a_protare.findInAncestry( summand->href( ) );
1498 if( ancestry == nullptr ) throw std::runtime_error( "Could not find href for summand - 1." );
1499 ancestry = ancestry->ancestor( );
1500 if( ancestry == nullptr ) throw std::runtime_error( "Could not find href for summand - 2." );
1501
1502 GIDI::Reaction const *GIDI_reaction2 = static_cast<GIDI::Reaction const *>( ancestry );
1503 for( std::size_t reactionIndex = 0; reactionIndex < m_reactions.size( ); ++reactionIndex ) {
1504 std::string label( m_reactions[reactionIndex]->label( ).c_str( ) );
1505
1506 if( label == GIDI_reaction2->label( ) ) {
1507 associatedOrphanProductIndices[reactionIndex].push_back( m_orphanProducts.size( ) - 1 );
1508 break;
1509 }
1510 }
1511 }
1512 }
1513
1514 for( std::size_t reactionIndex = 0; reactionIndex < m_reactions.size( ); ++reactionIndex ) {
1515 Reaction *reaction = m_reactions[reactionIndex];
1516 std::size_t size = associatedOrphanProductIndices[reactionIndex].size( );
1517 if( size > 0 ) {
1518 std::vector<Product *> associatedOrphanProducts;
1519 for( std::size_t index1 = 0; index1 < size; ++index1 ) {
1520 std::size_t associatedOrphanProductIndex = associatedOrphanProductIndices[reactionIndex][index1];
1521 m_orphanProducts[associatedOrphanProductIndex]->addOrphanProductToProductList( associatedOrphanProducts );
1522 }
1523 reaction->setOrphanProductData( associatedOrphanProductIndices[reactionIndex], associatedOrphanProducts );
1524 }
1525 }
1526 }
1527
1528 std::vector<GIDI::Reaction const *> GIDI_orphanProducts;
1529 for( std::size_t reactionIndex = 0; reactionIndex < a_protare.orphanProducts( ).size( ); ++reactionIndex ) {
1530 GIDI::Reaction const *GIDI_reaction = a_protare.orphanProduct( reactionIndex );
1531
1532 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue;
1533 GIDI_orphanProducts.push_back( GIDI_reaction );
1534 }
1535
1536 bool removeContinuousEnergyData = false;
1537 if( m_continuousEnergy ) {
1538 m_heatedCrossSections.update( a_smr, setupInfo, a_settings, particles, a_domainHash, a_temperatureInfos, GIDI_reactions, GIDI_orphanProducts,
1539 m_fixedGrid, zeroReactions );
1540 m_hasURR_probabilityTables = m_heatedCrossSections.hasURR_probabilityTables( );
1541 m_URR_domainMin = m_heatedCrossSections.URR_domainMin( );
1542 m_URR_domainMax = m_heatedCrossSections.URR_domainMax( ); }
1543 else {
1544 m_heatedMultigroupCrossSections.update( a_smr, a_protare, setupInfo, a_settings, particles, a_temperatureInfos, GIDI_reactions,
1545 GIDI_orphanProducts, zeroReactions, a_reactionsToExclude );
1546
1547 if( ( a_settings.upscatterModelAGroupBoundaries().size( ) > 0 ) ||
1548 ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::DBRC ) ) { // Load pointwise data to recompute Model A cross sections on user-defined grid or for DBRC.
1549 removeContinuousEnergyData = true;
1550 m_heatedCrossSections.update( a_smr, setupInfo, a_settings, particles, a_domainHash, a_temperatureInfos, GIDI_reactions, GIDI_orphanProducts,
1551 m_fixedGrid, zeroReactions );
1552 }
1553 }
1554
1555 if( m_upscatterModelASupported && ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::A ) ) {
1556 std::vector<double> const &upscatterModelAGroupBoundaries = a_settings.upscatterModelAGroupBoundaries( );
1557 if( upscatterModelAGroupBoundaries.size( ) == 0 ) {
1558 GIDI::Styles::Base const *style = a_protare.styles( ).get<GIDI::Styles::Base>( a_temperatureInfos[0].heatedMultiGroup( ) );
1559
1560 if( style->moniker( ) == GIDI_SnElasticUpScatterStyleChars ) style = a_protare.styles( ).get<GIDI::Styles::Base>( style->derivedStyle( ) );
1561 if( style->moniker( ) != GIDI_heatedMultiGroupStyleChars ) throw GIDI::Exception( "Label does not yield a heatedMultiGroup style." );
1562
1563 GIDI::Styles::HeatedMultiGroup const &heatedMultiGroup = *static_cast<GIDI::Styles::HeatedMultiGroup const *>( style );
1564 std::vector<double> const &boundaries = heatedMultiGroup.groupBoundaries( a_protare.projectile( ).ID( ) );
1565
1566 m_upscatterModelAGroupEnergies.resize( boundaries.size( ) );
1567 m_upscatterModelAGroupVelocities.resize( boundaries.size( ) );
1568 for( std::size_t i1 = 0; i1 < boundaries.size( ); ++i1 ) {
1569 m_upscatterModelAGroupEnergies[i1] = boundaries[i1];
1570 m_upscatterModelAGroupVelocities[i1] = MCGIDI_particleBeta( projectileMass( ), boundaries[i1] );
1571 }
1572
1573 GIDI::ExcludeReactionsSet reactionsToExclude;
1574 auto upscatterModelACrossSectionForm = a_protare.multiGroupCrossSection( a_smr, multiGroupSettings, a_temperatureInfos[0],
1575 reactionsToExclude, a_temperatureInfos[0].heatedMultiGroup( ) );
1576 m_upscatterModelACrossSection.resize( upscatterModelACrossSectionForm.size( ) );
1577 for( std::size_t i1 = 0; i1 < upscatterModelACrossSectionForm.size( ); ++i1 )
1578 m_upscatterModelACrossSection[i1] = upscatterModelACrossSectionForm[i1]; }
1579 else {
1580
1581 m_upscatterModelAGroupEnergies.reserve( upscatterModelAGroupBoundaries.size( ) );
1582 m_upscatterModelAGroupVelocities.reserve( upscatterModelAGroupBoundaries.size( ) );
1583 for( auto iter = upscatterModelAGroupBoundaries.begin( ); iter != upscatterModelAGroupBoundaries.end( ); ++iter ) {
1584 m_upscatterModelAGroupEnergies.push_back( *iter );
1585 m_upscatterModelAGroupVelocities.push_back( MCGIDI_particleBeta( projectileMass( ), *iter ) );
1586 }
1587
1588 GIDI::Transporting::MultiGroup boundaries( "Model A", upscatterModelAGroupBoundaries );
1589
1590 GIDI::Functions::XYs1d crossSectionXYs1d = m_heatedCrossSections.crossSectionAsGIDI_XYs1d( 0.0 );
1591
1592 GIDI::Transporting::Flux flux( "Model A", 0.0 );
1593 std::vector<double> energies, fluxes;
1594 energies.push_back( upscatterModelAGroupBoundaries[0] );
1595 energies.push_back( upscatterModelAGroupBoundaries.back( ) );
1596 fluxes.push_back( 1.0 );
1597 fluxes.push_back( 1.0 );
1598 flux.addFluxOrder( GIDI::Transporting::Flux_order( 0, energies, fluxes ) );
1599
1600 GIDI::Vector crossSectionVector = multiGroupXYs1d( boundaries, crossSectionXYs1d, flux );
1601 m_upscatterModelACrossSection.resize( crossSectionVector.size( ) );
1602 for( std::size_t index = 0; index < crossSectionVector.size( ); ++index )
1603 m_upscatterModelACrossSection[index] = crossSectionVector[index];
1604 }
1605 if( !m_continuousEnergy ) m_multiGroupHash = MultiGroupHash( m_projectileMultiGroupBoundariesCollapsed );
1606 }
1607
1608 if( ( PoPI::Intids::neutron == projectileIntid( ) ) && ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::DBRC ) ) {
1609 std::size_t reactionIndex = 0;
1610 for( auto reactionIter = m_reactions.begin( ); reactionIter != m_reactions.end( ); ++reactionIter, ++reactionIndex ) {
1611 if( (*reactionIter)->ENDF_MT( ) == 2 ) {
1612 Reaction *reaction = *reactionIter;
1613
1614 HeatedCrossSectionContinuousEnergy const *heatedCrossSectionContinuousEnergy = m_heatedCrossSections.heatedCrossSections( )[0];
1615 HeatedReactionCrossSectionContinuousEnergy const *heatedReactionCrossSectionContinuousEnergy
1616 = heatedCrossSectionContinuousEnergy->reactionCrossSection( reactionIndex );
1617
1618 Vector<double> const &energies = heatedCrossSectionContinuousEnergy->energies( );
1619 Vector<MCGIDI_FLOAT> const &crossSectionsFloat = heatedReactionCrossSectionContinuousEnergy->crossSections( );
1620 Vector<double> crossSections( crossSectionsFloat.size( ) );
1621 std::size_t index = 0;
1622 for( auto iter = crossSectionsFloat.begin( ); iter != crossSectionsFloat.end( ); ++iter, ++index )
1623 crossSections[index] = *iter;
1624
1625 Sampling::Upscatter::ModelDBRC_data *modelDBRC_data =
1626 new Sampling::Upscatter::ModelDBRC_data( projectileMass( ), targetMass( ), energies, crossSections, a_domainHash );
1627 reaction->setModelDBRC_data( modelDBRC_data );
1628 break;
1629 }
1630 }
1631 }
1632 if( removeContinuousEnergyData ) {
1633 m_heatedCrossSections.clear( );
1634 }
1635}
#define GIDI_multiGroupStyleChars
Definition GIDI.hpp:246
#define GIDI_MapInteractionAtomicChars
Definition GIDI.hpp:5164
#define GIDI_heatedMultiGroupStyleChars
Definition GIDI.hpp:253
#define GIDI_crossSectionSumChars
Definition GIDI.hpp:210
#define GIDI_SnElasticUpScatterStyleChars
Definition GIDI.hpp:254
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71
std::string const & label() const
Definition GIDI.hpp:658
bool active() const
Definition GIDI.hpp:4280
bool RutherfordScatteringPresent() const
Definition GIDI.hpp:4290
Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_label="") const
int ENDF_MT() const
Definition GIDI.hpp:4284
Component & crossSection()
Definition GIDI.hpp:4299
void productIDs(std::set< std::string > &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly) const
double crossSectionThreshold() const
Definition GIDI.hpp:4324
std::vector< double > groupBoundaries(std::string const &a_ID) const
GIDI::Suite const & transportables() const
Definition GIDI.hpp:3323
T * get(std::size_t a_Index)
Definition GIDI.hpp:2642
std::string const & href() const
Definition GIDI.hpp:4393
std::vector< double > const & boundaries() const
Definition GIDI.hpp:3501
MultiGroup multiGroup() const
Definition GIDI.hpp:3668
bool add(Particle const &a_particle)
bool hasParticle(std::string const &a_id) const
std::size_t size() const
Definition GIDI_data.hpp:79
Ancestry * findInAncestry(std::string const &a_href)
Ancestry * ancestor()
Definition GUPI.hpp:104
std::string const & moniker() const
Definition GUPI.hpp:102
statusMessageReporting * smr()
Definition LUPI.hpp:112
LUPI_HOST_DEVICE Reaction const * reaction(std::size_t a_index) const
Definition MCGIDI.hpp:1695
LUPI_HOST_DEVICE void crossSectionVector(double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const
LUPI_HOST_DEVICE int projectileIntid() const
Definition MCGIDI.hpp:1517
LUPI_HOST_DEVICE bool isTNSL_ProtareSingle() const
Definition MCGIDI.hpp:1541
LUPI_HOST_DEVICE double projectileMass() const
Definition MCGIDI.hpp:1520
LUPI_HOST_DEVICE double targetMass() const
Definition MCGIDI.hpp:1527
Vector multiGroupXYs1d(Transporting::MultiGroup const &a_boundaries, Functions::XYs1d const &a_function, Transporting::Flux const &a_flux)
Vector collapse(Vector const &a_vector, Transporting::Settings const &a_settings, Transporting::Particles const &a_particles, double a_temperature)
std::set< std::size_t > ExcludeReactionsSet
Definition GIDI.hpp:47
LUPI_HOST int MCGIDI_popsIntid(PoPI::Database const &a_pops, std::string const &a_ID)
LUPI_HOST int MCGIDI_popsIndex(PoPI::Database const &a_pops, std::string const &a_ID)
LUPI_HOST void convertACE_URR_probabilityTablesFromGIDI(GIDI::ProtareSingle const &a_protare, Transporting::MC &a_settings, SetupInfo &a_setupInfo)
#define smr_setReportError2p(smr, libraryID, code, fmt)
static std::string const photon
Definition PoPI.hpp:162
static std::string const electron
Definition PoPI.hpp:163
static int constexpr photon
Definition PoPI.hpp:178
static int constexpr neutron
Definition PoPI.hpp:177
static int constexpr electron
Definition PoPI.hpp:179

◆ ~ProtareSingle()

LUPI_HOST_DEVICE MCGIDI::ProtareSingle::~ProtareSingle ( )

Definition at line 1640 of file MCGIDI_protare.cc.

1640 {
1641
1642 for( Vector<NuclideGammaBranchInfo *>::const_iterator iter = m_branches.begin( ); iter < m_branches.end( ); ++iter ) delete *iter;
1643 for( Vector<NuclideGammaBranchStateInfo *>::const_iterator iter = m_nuclideGammaBranchStateInfos.begin( ); iter < m_nuclideGammaBranchStateInfos.end( ); ++iter ) delete *iter;
1644 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) delete *iter;
1645 for( Vector<Reaction *>::const_iterator iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) delete *iter;
1646}

Referenced by ~ProtareSingle().

Member Function Documentation

◆ branches()

LUPI_HOST_DEVICE const Vector< NuclideGammaBranchInfo * > & MCGIDI::ProtareSingle::branches ( ) const
inline

Returns a reference to the m_branches member.

Definition at line 1660 of file MCGIDI.hpp.

◆ continuousEnergy()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::continuousEnergy ( ) const
inline

Definition at line 1651 of file MCGIDI.hpp.

1651{ return( m_continuousEnergy ); }

◆ crossSection()

LUPI_HOST_DEVICE double MCGIDI::ProtareSingle::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 for target temperature a_temperature and projectile energy a_energy. a_sampling is only used for multi-group cross section look up.

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] The energy of the projectile.
a_sampling[in] Used for multi-group look up. If true, use augmented cross sections.

Definition at line 1818 of file MCGIDI_protare.cc.

1818 {
1819
1820 if( m_continuousEnergy ) return( m_heatedCrossSections.crossSection( a_URR_protareInfos, m_URR_index, a_hashIndex, a_temperature, a_energy ) );
1821
1822 return( m_heatedMultigroupCrossSections.crossSection( a_hashIndex, a_temperature, a_sampling ) );
1823}

Referenced by crossSection(), and sampleReaction().

◆ crossSectionVector()

LUPI_HOST_DEVICE void MCGIDI::ProtareSingle::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 1834 of file MCGIDI_protare.cc.

1835 {
1836
1837 if( m_continuousEnergy ) {
1838 if( !m_fixedGrid ) LUPI_THROW( "ProtareSingle::crossSectionVector: continuous energy cannot be supported." );
1839 m_heatedCrossSections.crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector ); }
1840 else {
1841 m_heatedMultigroupCrossSections.crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
1842 }
1843}
#define LUPI_THROW(arg)

Referenced by crossSectionVector().

◆ depositionEnergy()

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

Returns the index of a sampled reaction for a target with termpature a_temperature, a projectile with energy a_energy and total cross section a_crossSection. Random numbers are obtained via a_userrng and a_rngState.

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

Definition at line 1891 of file MCGIDI_protare.cc.

1891 {
1892
1893 if( m_continuousEnergy ) return( m_heatedCrossSections.depositionEnergy( a_hashIndex, a_temperature, a_energy ) );
1894
1895 return( m_heatedMultigroupCrossSections.depositionEnergy( a_hashIndex, a_temperature ) );
1896}

Referenced by depositionEnergy().

◆ depositionMomentum()

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

Returns the index of a sampled reaction for a target with termpature a_temperature, a projectile with energy a_energy and total cross section a_crossSection. Random numbers are obtained via a_userrng and a_rngState.

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

Definition at line 1907 of file MCGIDI_protare.cc.

1907 {
1908
1909 if( m_continuousEnergy ) return( m_heatedCrossSections.depositionMomentum( a_hashIndex, a_temperature, a_energy ) );
1910
1911 return( m_heatedMultigroupCrossSections.depositionMomentum( a_hashIndex, a_temperature ) );
1912}

Referenced by depositionMomentum().

◆ fixedGrid()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::fixedGrid ( ) const
inline

Definition at line 1652 of file MCGIDI.hpp.

1652{ return( m_fixedGrid ); }

◆ gain()

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

Returns the index of a sampled reaction for a target with termpature a_temperature, a projectile with energy a_energy and total cross section a_crossSection. Random numbers are obtained via a_userrng and a_rngState.

Parameters
a_hashIndex[in] Specifies the continuous energy or multi-group index.
a_temperature[in] The temperature of the target.
a_energy[in] The energy of the projectile.
a_particleIndex[in] The index of the particle whose gain is to be returned.

Definition at line 1940 of file MCGIDI_protare.cc.

1940 {
1941
1942 if( m_continuousEnergy ) return( m_heatedCrossSections.gain( a_hashIndex, a_temperature, a_energy, a_particleIndex ) );
1943
1944 return( m_heatedMultigroupCrossSections.gain( a_hashIndex, a_temperature, a_particleIndex ) );
1945}

Referenced by gain().

◆ gainViaIntid()

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

Returns the intid of a sampled reaction for a target with termpature a_temperature, a projectile with energy a_energy and total cross section a_crossSection. Random numbers are obtained via a_userrng and a_rngState.

Parameters
a_hashIndex[in] Specifies the continuous energy or multi-group index.
a_temperature[in] The temperature of the target.
a_energy[in] The energy of the projectile.
a_particleIntid[in] The intid of the particle whose gain is to be returned.

Definition at line 1957 of file MCGIDI_protare.cc.

1957 {
1958
1959 if( m_continuousEnergy ) return( m_heatedCrossSections.gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid ) );
1960
1961 return( m_heatedMultigroupCrossSections.gainViaIntid( a_hashIndex, a_temperature, a_particleIntid ) );
1962}

Referenced by gainViaIntid().

◆ hasFission()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::hasFission ( ) const

Returns true if this has a fission reaction and false otherwise.

Returns
true is if this has a fission reaction and false otherwise.

Definition at line 1771 of file MCGIDI_protare.cc.

1771 {
1772
1773 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) {
1774 if( (*iter)->hasFission( ) ) return( true );
1775 }
1776 return( false );
1777}

Referenced by hasFission().

◆ hasIncoherentDoppler()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::hasIncoherentDoppler ( ) const

Returns true if this has an incoherent photoatomic doppler broadened reaction and false otherwise.

Returns
true is if this has a specified reaction and false otherwise.

Definition at line 1785 of file MCGIDI_protare.cc.

1785 {
1786
1787 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) {
1788 if( (*iter)->ENDF_MT( ) == 1534 ) return( true );
1789 }
1790 return( false );
1791}

Referenced by hasIncoherentDoppler().

◆ hasURR_probabilityTables()

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

Definition at line 1706 of file MCGIDI.hpp.

1706{ return( m_hasURR_probabilityTables ); }

Referenced by MCGIDI::URR_protareInfos::setup().

◆ heatedCrossSections() [1/2]

LUPI_HOST_DEVICE HeatedCrossSectionsContinuousEnergy & MCGIDI::ProtareSingle::heatedCrossSections ( )
inline

Returns a reference to the m_heatedCrossSections member.

Definition at line 1654 of file MCGIDI.hpp.

◆ heatedCrossSections() [2/2]

LUPI_HOST_DEVICE HeatedCrossSectionsContinuousEnergy const & MCGIDI::ProtareSingle::heatedCrossSections ( ) const
inline

Returns a reference to the m_heatedCrossSections member.

Definition at line 1653 of file MCGIDI.hpp.

◆ heatedMultigroupCrossSections() [1/2]

LUPI_HOST_DEVICE HeatedCrossSectionsMultiGroup & MCGIDI::ProtareSingle::heatedMultigroupCrossSections ( )
inline

Returns a reference to the m_heatedMultigroupCrossSections member.

Definition at line 1656 of file MCGIDI.hpp.

◆ heatedMultigroupCrossSections() [2/2]

LUPI_HOST_DEVICE HeatedCrossSectionsMultiGroup const & MCGIDI::ProtareSingle::heatedMultigroupCrossSections ( ) const
inline

Returns a reference to the m_heatedMultigroupCrossSections member.

Definition at line 1655 of file MCGIDI.hpp.

◆ interaction()

LUPI_HOST_DEVICE String MCGIDI::ProtareSingle::interaction ( ) const
inline

Definition at line 1700 of file MCGIDI.hpp.

1700{ return( m_interaction ); }

Referenced by ProtareSingle().

◆ inURR()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::inURR ( double a_energy) const

Returns true if a_energy with unresolved resonance region (URR) of this and false otherwise.

Returns
true is if this has a URR data.

Definition at line 1799 of file MCGIDI_protare.cc.

1799 {
1800
1801 if( a_energy < m_URR_domainMin ) return( false );
1802 if( a_energy > m_URR_domainMax ) return( false );
1803
1804 return( true );
1805}

Referenced by inURR(), and MCGIDI::URR_protareInfos::updateProtare().

◆ isPhotoAtomic()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::isPhotoAtomic ( ) const
inline

Definition at line 1650 of file MCGIDI.hpp.

1650{ return( m_isPhotoAtomic ); }

Referenced by MCGIDI::HeatedCrossSectionContinuousEnergy::HeatedCrossSectionContinuousEnergy(), and ProtareSingle().

◆ maximumEnergy()

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

Returns the maximum cross section domain.

Definition at line 1684 of file MCGIDI.hpp.

◆ minimumEnergy()

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

Returns the minimum cross section domain.

Definition at line 1681 of file MCGIDI.hpp.

◆ nuclideGammaBranchStateInfos()

LUPI_HOST_DEVICE const Vector< NuclideGammaBranchStateInfo * > & MCGIDI::ProtareSingle::nuclideGammaBranchStateInfos ( ) const
inline

◆ numberOfOrphanProducts()

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

Returns the number of orphan products of this.

Definition at line 1696 of file MCGIDI.hpp.

◆ numberOfProtares()

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

Returns the number of protares contained in this.

Definition at line 1676 of file MCGIDI.hpp.

◆ numberOfReactions()

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

Returns the number of reactions of this.

Definition at line 1694 of file MCGIDI.hpp.

Referenced by protareWithReaction().

◆ orphanProduct()

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

Returns the (a_index-1)^th orphan product of this.

Definition at line 1697 of file MCGIDI.hpp.

◆ orphanProducts()

LUPI_HOST_DEVICE Vector< Reaction * > const & MCGIDI::ProtareSingle::orphanProducts ( ) const
inline

Returns the value of the m_orphanProducts member.

Definition at line 1666 of file MCGIDI.hpp.

◆ productionEnergy()

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

Returns the index of a sampled reaction for a target with termpature a_temperature, a projectile with energy a_energy and total cross section a_crossSection. Random numbers are obtained via a_userrng and a_rngState.

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

Definition at line 1923 of file MCGIDI_protare.cc.

1923 {
1924
1925 if( m_continuousEnergy ) return( m_heatedCrossSections.productionEnergy( a_hashIndex, a_temperature, a_energy ) );
1926
1927 return( m_heatedMultigroupCrossSections.productionEnergy( a_hashIndex, a_temperature ) );
1928}

Referenced by productionEnergy().

◆ projectileMultiGroupBoundaries()

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

Returns the value of the m_projectileMultiGroupBoundaries member.

Definition at line 1689 of file MCGIDI.hpp.

Referenced by MCGIDI::Protare::projectileMultiGroupBoundaries().

◆ projectileMultiGroupBoundariesCollapsed()

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

Returns the value of the m_projectileMultiGroupBoundariesCollapsed member.

Definition at line 1691 of file MCGIDI.hpp.

Referenced by MCGIDI::Protare::projectileMultiGroupBoundariesCollapsed(), and MCGIDI::HeatedCrossSectionsMultiGroup::update().

◆ protare() [1/2]

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

Returns the pointer representing the protare (i.e., this) if a_index is 0 and nullptr otherwise.

Parameters
a_index[in] Must always be 0.
Returns
Returns the pointer representing this.

Definition at line 1700 of file MCGIDI_protare.cc.

1700 {
1701
1702 if( a_index != 0 ) return( nullptr );
1703 return( this );
1704}

◆ protare() [2/2]

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

Returns the pointer representing the protare (i.e., this) if a_index is 0 and nullptr otherwise.

Parameters
a_index[in] Must always be 0.
Returns
Returns the pointer representing this.

Definition at line 1686 of file MCGIDI_protare.cc.

1686 {
1687
1688 if( a_index != 0 ) return( nullptr );
1689 return( this );
1690}

Referenced by protare(), and protare().

◆ protareWithReaction()

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

Returns the pointer to the this if (a_index - 1)th is a value reaction index and nullptr otherwise.

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

Definition at line 1714 of file MCGIDI_protare.cc.

1714 {
1715
1716 if( static_cast<std::size_t>( a_index ) < numberOfReactions( ) ) return( this );
1717 return( nullptr );
1718}
LUPI_HOST_DEVICE std::size_t numberOfReactions() const
Definition MCGIDI.hpp:1694

Referenced by protareWithReaction().

◆ reaction()

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

Returns the (a_index-1)^th reaction of this.

Definition at line 1695 of file MCGIDI.hpp.

Referenced by MCGIDI::HeatedCrossSectionContinuousEnergy::print().

◆ reactionCrossSection() [1/2]

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

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

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

Definition at line 1875 of file MCGIDI_protare.cc.

1875 {
1876
1877 if( m_continuousEnergy ) return( m_heatedCrossSections.reactionCrossSection( a_reactionIndex, a_URR_protareInfos, m_URR_index, a_temperature, a_energy ) );
1878
1879 return( m_heatedMultigroupCrossSections.reactionCrossSection( a_reactionIndex, a_temperature, a_energy ) );
1880}

◆ reactionCrossSection() [2/2]

LUPI_HOST_DEVICE double MCGIDI::ProtareSingle::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 reaction's cross section for the reaction at index a_reactionIndex, for target temperature a_temperature and projectile energy a_energy. a_sampling is only used for multi-group cross section look up.

Parameters
a_reactionIndex[in] The index of the reaction.
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] The energy of the projectile.
a_sampling[in] Used for multi-group look up. If true, use augmented cross sections.

Definition at line 1857 of file MCGIDI_protare.cc.

1858 {
1859
1860 if( m_continuousEnergy ) return( m_heatedCrossSections.reactionCrossSection( a_reactionIndex, a_URR_protareInfos, m_URR_index, a_hashIndex,
1861 a_temperature, a_energy ) );
1862
1863 return( m_heatedMultigroupCrossSections.reactionCrossSection( a_reactionIndex, a_hashIndex, a_temperature, a_sampling ) );
1864}

Referenced by MCGIDI::Distributions::CoherentPhotoAtomicScattering::angleBiasing(), MCGIDI::Distributions::IncoherentBoundToFreePhotoAtomicScattering::angleBiasing(), MCGIDI::Distributions::IncoherentPhotoAtomicScattering::angleBiasing(), reactionCrossSection(), and reactionCrossSection().

◆ reactionHasURR_probabilityTables()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::reactionHasURR_probabilityTables ( std::size_t a_index) const
inline

Definition at line 1709 of file MCGIDI.hpp.

1709{ return( m_heatedCrossSections.reactionHasURR_probabilityTables( a_index ) ); }

◆ reactions()

LUPI_HOST_DEVICE Vector< Reaction * > const & MCGIDI::ProtareSingle::reactions ( ) const
inline

Returns the value of the m_reactions member.

Definition at line 1664 of file MCGIDI.hpp.

◆ sampleBranchingGammas()

template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void MCGIDI::ProtareSingle::sampleBranchingGammas ( Sampling::Input & a_input,
double a_projectileEnergy,
int a_initialStateIndex,
RNG && a_rng,
PUSHBACK && a_push_back,
Sampling::ProductHandler & a_products ) const
inline

Samples gammas from a nuclide electro-magnetic decay.

Parameters
a_input[in] Sample options requested by user.
a_projectileEnergy[in] The energy of the projectile.
a_initialStateIndex[in] The index in m_nuclideGammaBranchStateInfos whose nuclide data are used for sampling.
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 gammas to.

Definition at line 2718 of file MCGIDI_headerSource.hpp.

2719 {
2720
2721 int initialStateIndex = a_initialStateIndex;
2722 double energyLevelSampleWidthUpper = 0.0; // Used for GRIN continuum levels to add variaction to outgoing photons.
2723
2724 NuclideGammaBranchStateInfo *nuclideGammaBranchStateInfo = nullptr;
2725 if( initialStateIndex >= 0 ) nuclideGammaBranchStateInfo = m_nuclideGammaBranchStateInfos[static_cast<std::size_t>(initialStateIndex)];
2726 while( initialStateIndex >= 0 ) {
2727 auto const &branchIndices = nuclideGammaBranchStateInfo->branchIndices( );
2728
2729 double random = a_rng( );
2730 double sum = 0.0;
2731 initialStateIndex = -1; // Just in case the for loop never has "sum >= random".
2732 for( std::size_t i1 = 0; i1 < branchIndices.size( ); ++i1 ) {
2733 NuclideGammaBranchInfo *nuclideGammaBranchInfo = m_branches[branchIndices[i1]];
2734
2735 sum += nuclideGammaBranchInfo->probability( );
2736 if( sum >= random ) {
2737 double energyLevelSampleWidthLower = 0.0;
2738 initialStateIndex = nuclideGammaBranchInfo->residualStateIndex( );
2739 if( initialStateIndex >= 0 ) {
2740 nuclideGammaBranchStateInfo = m_nuclideGammaBranchStateInfos[static_cast<std::size_t>(initialStateIndex)];
2741 energyLevelSampleWidthLower = a_rng( ) * nuclideGammaBranchStateInfo->nuclearLevelEnergyWidth( );
2742 }
2743 if( nuclideGammaBranchInfo->photonEmissionProbability( ) > a_rng( ) ) {
2744 a_input.setSampledType( Sampling::SampledType::photon );
2745 a_input.m_dataInTargetFrame = false;
2746 a_input.m_frame = GIDI::Frame::lab;
2747
2748 a_input.m_energyOut1 = nuclideGammaBranchInfo->gammaEnergy( ) + energyLevelSampleWidthUpper - energyLevelSampleWidthLower;
2749 a_input.m_mu = 1.0 - 2.0 * a_rng( );
2750 a_input.m_phi = 2.0 * M_PI * a_rng( );
2751
2752 a_products.add( a_projectileEnergy, PoPI::Intids::photon, m_photonIndex, userPhotonIndex( ), 0.0, a_input, a_rng, a_push_back, true );
2753 }
2754 energyLevelSampleWidthUpper = energyLevelSampleWidthLower;
2755 break;
2756 }
2757 }
2758 }
2759}
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE int userPhotonIndex() const
Definition MCGIDI.hpp:1531

Referenced by MCGIDI::Product::sampleFinalState(), MCGIDI::GRIN_capture::sampleProducts(), MCGIDI::GRIN_inelastic::sampleProducts(), and MCGIDI::Product::sampleProducts().

◆ sampleReaction()

template<typename RNG>
LUPI_HOST_DEVICE std::size_t MCGIDI::ProtareSingle::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 index of a sampled reaction for target temperature, projectile energy and total cross section as specified via argument a_input. Random numbers are obtained via a_rng.

Parameters
a_input[in/out] Sample options requested by user. The values m_modelTemperature and m_modelEnergy are set by this method.
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).

Definition at line 2843 of file MCGIDI_headerSource.hpp.

2844 {
2845
2846 std::size_t hashIndex = a_hashIndex;
2847 double crossSection1 = a_crossSection;
2848
2849 a_input.m_dataInTargetFrame = false;
2850 a_input.m_modelTemperature = a_input.m_temperature;
2851 a_input.m_modelEnergy = a_input.m_energy;
2852
2853 if( upscatterModelASupported( ) && ( a_input.m_upscatterModel == Sampling::Upscatter::Model::A ) ) {
2854 a_input.m_dataInTargetFrame = sampleTargetBetaForUpscatterModelA( a_input, a_rng );
2855 if( a_input.m_dataInTargetFrame ) {
2856 if( m_continuousEnergy ) {
2857 hashIndex = m_domainHash.index( a_input.m_modelEnergy ); }
2858 else {
2859 hashIndex = m_multiGroupHash.index( a_input.m_modelEnergy );
2860 }
2861 crossSection1 = crossSection( a_URR_protareInfos, hashIndex, a_input.m_modelTemperature, a_input.m_modelEnergy, true );
2862 }
2863 }
2864
2865 if( m_continuousEnergy ) return( m_heatedCrossSections.sampleReaction( a_URR_protareInfos, m_URR_index, hashIndex,
2866 a_input.m_modelTemperature, a_input.m_modelEnergy, crossSection1, a_rng ) );
2867
2868 return( m_heatedMultigroupCrossSections.sampleReaction( hashIndex, a_input.m_modelTemperature, a_input.m_modelEnergy,
2869 crossSection1, a_rng ) );
2870}
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 bool sampleTargetBetaForUpscatterModelA(Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE bool upscatterModelASupported() const
Definition MCGIDI.hpp:1735

◆ sampleTargetBetaForUpscatterModelA()

template<typename RNG>
LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::sampleTargetBetaForUpscatterModelA ( Sampling::Input & a_input,
RNG && a_rng ) const
inline

This function is used internally to sample a target's velocity (speed and cosine of angle relative to projectile) for a heated target using zero temperature, multi-grouped cross sections.

Parameters
a_input[in/out] Contains needed input like the targets temperature. Also will have the target sampled velocity on return if return value is true.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
Returns
Returns true if target velocity is sampled and false otherwise.

Definition at line 2772 of file MCGIDI_headerSource.hpp.

2772 {
2773
2774 double projectileBeta = MCGIDI_particleBeta( m_projectileMass, a_input.energy( ) );
2775 double targetThermalBeta = MCGIDI_particleBeta( m_targetMass, a_input.temperature( ) );
2776
2777 a_input.m_projectileBeta = projectileBeta;
2778 a_input.m_relativeBeta = projectileBeta;
2779 a_input.m_muLab = 0.0;
2780 a_input.m_targetBeta = 0.0;
2781
2782 if( targetThermalBeta < 1e-4 * projectileBeta ) return( false );
2783
2784 a_input.m_modelTemperature = 0.0;
2785
2786 double relativeBetaMin = projectileBeta - 2.0 * targetThermalBeta;
2787 double relativeBetaMax = projectileBeta + 2.0 * targetThermalBeta;
2788
2789 std::size_t maxIndex = m_upscatterModelAGroupVelocities.size( ) - 2;
2790 int intRelativeBetaMinIndex = binarySearchVector( relativeBetaMin, m_upscatterModelAGroupVelocities, true );
2791 std::size_t relativeBetaMinIndex = static_cast<std::size_t>( intRelativeBetaMinIndex );
2792 int intRelativeBetaMaxIndex = binarySearchVector( relativeBetaMax, m_upscatterModelAGroupVelocities, true );
2793 std::size_t relativeBetaMaxIndex = static_cast<std::size_t>( intRelativeBetaMaxIndex );
2794 double targetBeta, relativeBeta, mu;
2795
2796 if( relativeBetaMinIndex >= maxIndex ) relativeBetaMinIndex = maxIndex;
2797 if( relativeBetaMaxIndex >= maxIndex ) relativeBetaMaxIndex = maxIndex;
2798
2799 if( relativeBetaMinIndex == relativeBetaMaxIndex ) {
2800 targetBeta = targetThermalBeta * sampleBetaFromMaxwellian( a_rng );
2801 mu = 1.0 - 2.0 * a_rng( );
2802 relativeBeta = sqrt( targetBeta * targetBeta + projectileBeta * projectileBeta - 2.0 * mu * targetBeta * projectileBeta ); }
2803 else {
2804
2805 double reactionRate;
2806 double reactionRateMax = 0;
2807 for( std::size_t i1 = relativeBetaMinIndex; i1 <= relativeBetaMaxIndex; ++i1 ) {
2808 reactionRate = m_upscatterModelACrossSection[i1] * m_upscatterModelAGroupVelocities[i1+1];
2809 if( reactionRate > reactionRateMax ) reactionRateMax = reactionRate;
2810 }
2811
2812 do {
2813 targetBeta = targetThermalBeta * sampleBetaFromMaxwellian( a_rng );
2814 mu = 1.0 - 2.0 * a_rng( );
2815 relativeBeta = sqrt( targetBeta * targetBeta + projectileBeta * projectileBeta - 2.0 * mu * targetBeta * projectileBeta );
2816
2817 std::size_t index = static_cast<std::size_t>( binarySearchVector( relativeBeta, m_upscatterModelAGroupVelocities, true ) );
2818 if( index > maxIndex ) index = maxIndex;
2819 reactionRate = m_upscatterModelACrossSection[index] * relativeBeta;
2820 } while( reactionRate < a_rng( ) * reactionRateMax );
2821 }
2822
2823 a_input.m_modelEnergy = particleKineticEnergy( m_projectileMass, relativeBeta );
2824 a_input.m_relativeBeta = relativeBeta;
2825 a_input.m_muLab = mu;
2826 a_input.m_targetBeta = targetBeta;
2827
2828 return( true );
2829}
LUPI_HOST_DEVICE double sampleBetaFromMaxwellian(RNG &&a_rng)
LUPI_HOST_DEVICE double particleKineticEnergy(double a_mass_unitOfEnergy, double a_particleBeta)
LUPI_HOST_DEVICE int binarySearchVector(double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
Definition MCGIDI.hpp:318

Referenced by sampleReaction().

◆ serialize2()

LUPI_HOST_DEVICE void MCGIDI::ProtareSingle::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 1972 of file MCGIDI_protare.cc.

1972 {
1973
1974 std::size_t vectorSize;
1975 LUPI::DataBuffer *workingBuffer = &a_buffer;
1976
1977 DATA_MEMBER_STRING( m_interaction, a_buffer, a_mode );
1978 DATA_MEMBER_INT( m_URR_index, a_buffer, a_mode );
1979 DATA_MEMBER_CAST( m_hasURR_probabilityTables, a_buffer, a_mode, bool );
1980 DATA_MEMBER_DOUBLE( m_URR_domainMin, a_buffer, a_mode );
1981 DATA_MEMBER_DOUBLE( m_URR_domainMax, a_buffer, a_mode );
1982 m_domainHash.serialize( a_buffer, a_mode );
1983 DATA_MEMBER_VECTOR_DOUBLE( m_projectileMultiGroupBoundaries, a_buffer, a_mode );
1984 DATA_MEMBER_VECTOR_DOUBLE( m_projectileMultiGroupBoundariesCollapsed, a_buffer, a_mode );
1985 DATA_MEMBER_CAST( m_upscatterModelASupported, a_buffer, a_mode, bool );
1986 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelAGroupEnergies, a_buffer, a_mode );
1987 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelAGroupVelocities, a_buffer, a_mode );
1988 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelACrossSection, a_buffer, a_mode );
1989 m_multiGroupHash.serialize( a_buffer, a_mode );
1990
1991 vectorSize = m_nuclideGammaBranchStateInfos.size( );
1992 int vectorSizeInt = (int) vectorSize;
1993 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
1994 vectorSize = (std::size_t) vectorSizeInt;
1995
1996 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
1997 m_nuclideGammaBranchStateInfos.resize( vectorSize, &(workingBuffer->m_placement) );
1998 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
1999 if (workingBuffer->m_placement != nullptr) {
2000 m_nuclideGammaBranchStateInfos[vectorIndex] = new(workingBuffer->m_placement) NuclideGammaBranchStateInfo;
2001 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchStateInfo ) );
2002 }
2003 else {
2004 m_nuclideGammaBranchStateInfos[vectorIndex] = new NuclideGammaBranchStateInfo;
2005 }
2006 }
2007 }
2008 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2009 a_buffer.m_placement += m_nuclideGammaBranchStateInfos.internalSize();
2010 a_buffer.incrementPlacement( sizeof( NuclideGammaBranchStateInfo ) * vectorSize );
2011 }
2012 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2013 m_nuclideGammaBranchStateInfos[vectorIndex]->serialize( *workingBuffer, a_mode );
2014 }
2015
2016 vectorSize = m_branches.size( );
2017 vectorSizeInt = (int) vectorSize;
2018 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2019 vectorSize = (std::size_t) vectorSizeInt;
2020
2021 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2022 m_branches.resize( vectorSize, &(workingBuffer->m_placement) );
2023 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2024 if (workingBuffer->m_placement != nullptr) {
2025 m_branches[vectorIndex] = new(workingBuffer->m_placement) NuclideGammaBranchInfo;
2026 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchInfo ) );
2027 }
2028 else {
2029 m_branches[vectorIndex] = new NuclideGammaBranchInfo;
2030 }
2031 }
2032 }
2033 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2034 a_buffer.m_placement += m_branches.internalSize();
2035 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchInfo ) * vectorSize );
2036 }
2037 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2038 m_branches[vectorIndex]->serialize( *workingBuffer, a_mode );
2039 }
2040
2041 vectorSize = m_reactions.size( );
2042 vectorSizeInt = (int) vectorSize;
2043 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2044 vectorSize = (std::size_t) vectorSizeInt;
2045
2046 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2047 m_reactions.resize( vectorSize, &(workingBuffer->m_placement) );
2048 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2049 if (workingBuffer->m_placement != nullptr) {
2050 m_reactions[vectorIndex] = new(workingBuffer->m_placement) Reaction;
2051 workingBuffer->incrementPlacement( sizeof(Reaction));
2052 }
2053 else {
2054 m_reactions[vectorIndex] = new Reaction;
2055 }
2056 }
2057 }
2058 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2059 a_buffer.m_placement += m_reactions.internalSize();
2060 a_buffer.incrementPlacement( sizeof(Reaction) * vectorSize);
2061 }
2062 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2063 m_reactions[vectorIndex]->serialize( *workingBuffer, a_mode );
2064 m_reactions[vectorIndex]->updateProtareSingleInfo( this, vectorIndex );
2065 }
2066
2067 vectorSize = m_orphanProducts.size( );
2068 vectorSizeInt = (int) vectorSize;
2069 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2070 vectorSize = (std::size_t) vectorSizeInt;
2071
2072 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2073 m_orphanProducts.resize( vectorSize, &(workingBuffer->m_placement) );
2074 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2075 if (workingBuffer->m_placement != nullptr) {
2076 m_orphanProducts[vectorIndex] = new(workingBuffer->m_placement) Reaction;
2077 workingBuffer->incrementPlacement( sizeof(Reaction));
2078 }
2079 else {
2080 m_orphanProducts[vectorIndex] = new Reaction;
2081 }
2082 }
2083 }
2084
2085 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2086 a_buffer.m_placement += m_orphanProducts.internalSize( );
2087 a_buffer.incrementPlacement( sizeof( Reaction ) * vectorSize );
2088 }
2089
2090 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2091 m_orphanProducts[vectorIndex]->serialize( *workingBuffer, a_mode );
2092 m_orphanProducts[vectorIndex]->updateProtareSingleInfo( this, vectorIndex );
2093 }
2094
2095 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2096 for( auto reactionIter = m_reactions.begin( ); reactionIter != m_reactions.end( ); ++reactionIter ) {
2097 (*reactionIter)->addOrphanProductToProductList( m_orphanProducts );
2098 }
2099 }
2100
2101 DATA_MEMBER_CAST( m_isPhotoAtomic, *workingBuffer, a_mode, bool );
2102 DATA_MEMBER_CAST( m_continuousEnergy, *workingBuffer, a_mode, bool );
2103 DATA_MEMBER_CAST( m_fixedGrid, *workingBuffer, a_mode, bool );
2104 m_heatedCrossSections.serialize( *workingBuffer, a_mode );
2105 m_heatedMultigroupCrossSections.serialize( *workingBuffer, a_mode );
2106}
#define DATA_MEMBER_STRING(member, buf, mode)
#define DATA_MEMBER_CAST(member, buf, mode, someType)
#define DATA_MEMBER_VECTOR_DOUBLE(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)

Referenced by serialize2().

◆ setURR_index()

LUPI_HOST_DEVICE void MCGIDI::ProtareSingle::setURR_index ( int a_URR_index)
inline

Definition at line 1704 of file MCGIDI.hpp.

1704{ m_URR_index = a_URR_index; }

Referenced by MCGIDI::URR_protareInfos::setup().

◆ setUserParticleIndex2()

LUPI_HOST void MCGIDI::ProtareSingle::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 1655 of file MCGIDI_protare.cc.

1655 {
1656
1657 m_heatedCrossSections.setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1658 m_heatedMultigroupCrossSections.setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1659 for( auto iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1660 for( auto iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1661}

Referenced by setUserParticleIndex2().

◆ setUserParticleIndexViaIntid2()

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

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

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 1670 of file MCGIDI_protare.cc.

1670 {
1671
1672 m_heatedCrossSections.setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1673 m_heatedMultigroupCrossSections.setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1674 for( auto iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1675 for( auto iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1676}

Referenced by setUserParticleIndexViaIntid2().

◆ sizeOf2()

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

Definition at line 1742 of file MCGIDI.hpp.

1742{ return sizeof(*this); }

◆ temperatures()

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

Returns the list of temperatures for this.

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

Definition at line 1728 of file MCGIDI_protare.cc.

1728 {
1729
1730 if( a_index != 0 ) LUPI_THROW( "ProtareSingle::temperatures: a_index not 0." );
1731 if( m_continuousEnergy ) return( m_heatedCrossSections.temperatures( ) );
1732 return( m_heatedMultigroupCrossSections.temperatures( ) );
1733}

Referenced by MCGIDI::Protare::temperatures(), and temperatures().

◆ threshold()

LUPI_HOST_DEVICE double MCGIDI::ProtareSingle::threshold ( std::size_t a_index) const
inline

Returns the threshold for the reaction at index a_index.

Definition at line 1711 of file MCGIDI.hpp.

◆ upscatterModelACrossSection()

LUPI_HOST_DEVICE Vector< double > const & MCGIDI::ProtareSingle::upscatterModelACrossSection ( ) const
inline

Returns the value of the m_upscatterModelACrossSection.

Definition at line 1738 of file MCGIDI.hpp.

◆ upscatterModelAGroupEnergies()

LUPI_HOST_DEVICE Vector< double > const & MCGIDI::ProtareSingle::upscatterModelAGroupEnergies ( ) const
inline

Returns a reference to the m_upscatterModelAGroupEnergies member.

Definition at line 1736 of file MCGIDI.hpp.

◆ upscatterModelAGroupVelocities()

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

Returns a reference to the m_upscatterModelAGroupVelocities member.

Definition at line 1737 of file MCGIDI.hpp.

◆ upscatterModelASupported()

LUPI_HOST_DEVICE bool MCGIDI::ProtareSingle::upscatterModelASupported ( ) const
inline

Returns the value of the m_upscatterModelASupported member.

Definition at line 1735 of file MCGIDI.hpp.

Referenced by sampleReaction().

◆ URR_domainMax()

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

Definition at line 1708 of file MCGIDI.hpp.

1708{ return( m_URR_domainMax ); }

◆ URR_domainMin()

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

Definition at line 1707 of file MCGIDI.hpp.

1707{ return( m_URR_domainMin ); }

◆ URR_index()

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

Definition at line 1703 of file MCGIDI.hpp.

1703{ return( m_URR_index ); }

Referenced by MCGIDI::URR_protareInfos::updateProtare().


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