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

#include <MCGIDI.hpp>

Public Member Functions

LUPI_HOST_DEVICE HeatedCrossSectionContinuousEnergy ()
LUPI_HOST HeatedCrossSectionContinuousEnergy (SetupInfo &a_setupInfo, Transporting::MC const &a_settings, GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash, GIDI::Styles::TemperatureInfo const &a_temperatureInfo, std::vector< GIDI::Reaction const * > const &a_reactions, std::vector< GIDI::Reaction const * > const &a_orphanProducts, bool a_fixedGrid, bool a_zeroReactions)
LUPI_HOST_DEVICE ~HeatedCrossSectionContinuousEnergy ()
LUPI_HOST_DEVICE std::size_t evaluationInfo (std::size_t a_hashIndex, double a_energy, double *a_energyFraction) const
LUPI_HOST HeatedReactionCrossSectionContinuousEnergy const * reactionCrossSection (std::size_t a_index) const
LUPI_HOST_DEVICE double temperature () const
LUPI_HOST_DEVICE double minimumEnergy () const
LUPI_HOST_DEVICE double maximumEnergy () const
LUPI_HOST_DEVICE std::size_t numberOfReactions () const
LUPI_HOST_DEVICE std::size_t thresholdOffset (std::size_t a_reactionIndex) const
LUPI_HOST_DEVICE double threshold (std::size_t a_reactionIndex) 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 Vector< MCGIDI_FLOAT > & totalCrossSection ()
LUPI_HOST_DEVICE double crossSection (URR_protareInfos const &a_URR_protareInfos, int a_URR_index, std::size_t a_hashIndex, double a_energy, bool a_sampling=false) const
LUPI_HOST GIDI::Functions::XYs1d crossSectionAsGIDI_XYs1d () const
LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, int a_URR_index, std::size_t a_hashIndex, double a_energy, bool a_sampling=false) const
LUPI_HOST_DEVICE double reactionCrossSection2 (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, int a_URR_index, double a_energy, std::size_t a_energyIndex, double a_energyFraction, bool a_sampling=false) const
LUPI_HOST_DEVICE double reactionCrossSection (std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, int a_URR_index, double a_energy) const
LUPI_HOST GIDI::Functions::XYs1d reactionCrossSectionAsGIDI_XYs1d (std::size_t a_reactionIndex) const
LUPI_HOST_DEVICE double depositionEnergy (std::size_t a_hashIndex, double a_energy) const
LUPI_HOST_DEVICE double depositionMomentum (std::size_t a_hashIndex, double a_energy) const
LUPI_HOST_DEVICE double productionEnergy (std::size_t a_hashIndex, double a_energy) const
LUPI_HOST_DEVICE double gain (std::size_t a_hashIndex, double a_energy, int a_particleIndex) const
LUPI_HOST_DEVICE double gainViaIntid (std::size_t a_hashIndex, double a_energy, int a_particleIntid) const
LUPI_HOST void setUserParticleIndex (int a_particleIndex, int a_userParticleIndex)
LUPI_HOST void setUserParticleIndexViaIntid (int a_particleIntid, int a_userParticleIndex)
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Vector< double > const & energies () const
LUPI_HOST void print (ProtareSingle const *a_protareSingle, std::string const &a_indent, std::string const &a_iFormat, std::string const &a_energyFormat, std::string const &a_dFormat) const

Detailed Description

Definition at line 614 of file MCGIDI.hpp.

Constructor & Destructor Documentation

◆ HeatedCrossSectionContinuousEnergy() [1/2]

LUPI_HOST_DEVICE MCGIDI::HeatedCrossSectionContinuousEnergy::HeatedCrossSectionContinuousEnergy ( )

Definition at line 268 of file MCGIDI_heatedCrossSections.cc.

268 :
269 m_temperature( 0.0 ),
270 m_hashIndices( ),
271 m_energies( ),
272 m_totalCrossSection( ),
273 m_depositionEnergy( ),
274 m_depositionMomentum( ),
275 m_productionEnergy( ),
276 m_gains( ),
277 m_URR_mode( Transporting::URR_mode::none ),
278 m_reactionsInURR_region( ),
279 m_reactionCrossSections( ),
280 m_ACE_URR_probabilityTables( nullptr ) {
281
282}

◆ HeatedCrossSectionContinuousEnergy() [2/2]

LUPI_HOST MCGIDI::HeatedCrossSectionContinuousEnergy::HeatedCrossSectionContinuousEnergy ( SetupInfo & a_setupInfo,
Transporting::MC const & a_settings,
GIDI::Transporting::Particles const & a_particles,
DomainHash const & a_domainHash,
GIDI::Styles::TemperatureInfo const & a_temperatureInfo,
std::vector< GIDI::Reaction const * > const & a_reactions,
std::vector< GIDI::Reaction const * > const & a_orphanProducts,
bool a_fixedGrid,
bool a_zeroReactions )

Fills in this with the requested temperature data.

Parameters
a_setupInfo[in] Used internally when constructing a Protare to pass information to other components.
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_temperatureInfo[in] The list of temperatures to use.
a_reactions[in] The list of reactions to use.
a_orphanProducts[in] The list of orphan products to use.
a_fixedGrid[in] If true, the specified fixed grid is used; otherwise, grid in the file is used.
a_zeroReactions[in] Special case where no reaction in a protare is wanted so the first one is used but its cross section is set to 0.0 at all energies.

Definition at line 298 of file MCGIDI_heatedCrossSections.cc.

301 :
302 m_temperature( a_temperatureInfo.temperature( ).value( ) ),
303 m_hashIndices( ),
304 m_energies( ),
305 m_totalCrossSection( ),
306 m_depositionEnergy( ),
307 m_depositionMomentum( ),
308 m_productionEnergy( ),
309 m_gains( ),
310 m_URR_mode( Transporting::URR_mode::none ),
311 m_reactionsInURR_region( ),
312 m_reactionCrossSections( ),
313 m_ACE_URR_probabilityTables( nullptr ) {
314
315 bool isPhotoAtomic = a_setupInfo.m_protare.isPhotoAtomic( );
316 std::string label( a_temperatureInfo.griddedCrossSection( ) );
317 std::string URR_label( a_temperatureInfo.URR_probabilityTables( ) );
318
319 GIDI::Styles::GriddedCrossSection const &griddedCrossSectionStyle =
320 static_cast<GIDI::Styles::GriddedCrossSection const &>( *a_settings.styles( )->get<GIDI::Styles::Base>( label ) );
321 GIDI::Grid const &grid = griddedCrossSectionStyle.grid( );
322
323 std::vector<double> const *energiesPointer;
324 std::vector<double> const &energies = grid.data( ).vector();
325 std::vector<double> const &fixedGridPoints = a_settings.fixedGridPoints( );
326 std::vector<int> fixedGridIndices( fixedGridPoints.size( ) );
327 if( a_fixedGrid ) {
328 for( std::size_t i1 = 0; i1 < fixedGridPoints.size( ); ++i1 ) {
329 fixedGridIndices[i1] = binarySearchVector( fixedGridPoints[i1], energies );
330 }
331 energiesPointer = &fixedGridPoints;
332 m_energies = fixedGridPoints; }
333 else {
334 energiesPointer = &energies;
335 m_energies = energies;
336 }
337
338 m_hashIndices = a_domainHash.map( m_energies );
339
340 std::size_t reactionIndex = 0;
341 GIDI::Axes axes;
342 std::vector<double> dummy;
343 GIDI::Functions::Ys1d totalCrossSection( axes, ptwXY_interpolationLinLin, 0, dummy );
344 GIDI::Functions::Ys1d fixedGridCrossSection( axes, ptwXY_interpolationLinLin, 0, fixedGridPoints );
345 m_reactionCrossSections.resize( a_reactions.size( ) );
346 m_URR_mode = Transporting::URR_mode::none;
347 for( std::vector<GIDI::Reaction const *>::const_iterator reactionIter = a_reactions.begin( ); reactionIter != a_reactions.end( ); ++reactionIter, ++reactionIndex ) {
348 GIDI::Suite const &reactionCrossSectionSuite = (*reactionIter)->crossSection( );
349 GIDI::Functions::Ys1d const *reactionCrossSection3 = reactionCrossSectionSuite.get<GIDI::Functions::Ys1d>( label );
350 GIDI::Functions::Ys1d *reactionCrossSectionZeroReactions = nullptr;
351 if( a_zeroReactions ) {
352 reactionCrossSectionZeroReactions = new GIDI::Functions::Ys1d( *reactionCrossSection3 );
353 for( std::size_t index = 0; index < reactionCrossSectionZeroReactions->size( ); ++index ) reactionCrossSectionZeroReactions->set( index, 0.0 );
354 reactionCrossSection3 = reactionCrossSectionZeroReactions;
355 }
356
357 Probabilities::ProbabilityBase2d *URR_probabilityTables = nullptr;
358 if( ( a_settings._URR_mode( ) == Transporting::URR_mode::pdfs ) && ( URR_label != "" ) ) {
359 if( reactionCrossSectionSuite.has( URR_label ) ) {
360 GIDI::Functions::URR_probabilityTables1d const &URR_probability_tables1d( *reactionCrossSectionSuite.get<GIDI::Functions::URR_probabilityTables1d>( URR_label ) );
361 URR_probabilityTables = Probabilities::parseProbability2d( URR_probability_tables1d.function2d( ), nullptr );
362 m_URR_mode = Transporting::URR_mode::pdfs;
363 }
364 }
365
366 ACE_URR_probabilityTables *ACE_URR_probabilityTables1 = nullptr;
367 if( a_settings._URR_mode( ) == Transporting::URR_mode::ACE_URR_probabilityTables ) {
368 auto URR_iter = a_setupInfo.m_ACE_URR_probabilityTablesFromGIDI.find( URR_label );
369 if( URR_iter != a_setupInfo.m_ACE_URR_probabilityTablesFromGIDI.end( ) ) {
370 auto ACE_URR_probabilityTablesIter = (*URR_iter).second->m_ACE_URR_probabilityTables.find( (*reactionIter)->label( ) );
371 if( ACE_URR_probabilityTablesIter != (*URR_iter).second->m_ACE_URR_probabilityTables.end( ) ) {
372 ACE_URR_probabilityTables1 = (*ACE_URR_probabilityTablesIter).second;
373 (*ACE_URR_probabilityTablesIter).second = nullptr; // Set to nullptr so destructor of m_ACE_URR_probabilityTables does not delete.
375 }
376 }
377 }
378
379 if( a_fixedGrid ) {
380 GIDI::Functions::Ys1d *reactionCrossSection4 = &fixedGridCrossSection;
381
382 std::size_t start = 0;
383
384 if( energies[reactionCrossSection3->start( )] > fixedGridPoints[0] ) {
385 int intStart = binarySearchVector( energies[reactionCrossSection3->start( )], fixedGridPoints ) + 1;
386 if (intStart < 0) {
387 throw std::out_of_range( "Cross section out of range when initializing fixed grid!" );
388 }
389 start = static_cast<std::size_t>( intStart );
390 }
391
392 for( std::size_t i1 = 0; i1 < start; ++i1 ) reactionCrossSection4->set( i1, 0.0 );
393 for( std::size_t i1 = start; i1 < fixedGridPoints.size( ); ++i1 ) {
394 std::size_t index = static_cast<std::size_t>( fixedGridIndices[i1] );
395 double fraction = ( fixedGridPoints[i1] - energies[index] ) / ( energies[index+1] - energies[index] );
396
397 index -= reactionCrossSection3->start( );
398 reactionCrossSection4->set( i1, ( 1.0 - fraction ) * (*reactionCrossSection3)[index] + fraction * (*reactionCrossSection3)[index+1] );
399 }
400 reactionCrossSection3 = &fixedGridCrossSection;
401 }
402
403 m_reactionCrossSections[reactionIndex] = new HeatedReactionCrossSectionContinuousEnergy( (*reactionIter)->crossSectionThreshold( ),
404 *reactionCrossSection3, URR_probabilityTables, ACE_URR_probabilityTables1 );
405 totalCrossSection += *reactionCrossSection3;
406
407 delete reactionCrossSectionZeroReactions;
408 }
409 m_totalCrossSection.resize( totalCrossSection.length( ), 0.0 );
410 for( std::size_t i1 = 0; i1 < totalCrossSection.size( ); ++i1 ) m_totalCrossSection[i1+totalCrossSection.start()] = totalCrossSection[i1];
411
412 if( hasURR_probabilityTables( ) ) {
413 std::vector<std::size_t> reactions_in_URR_region;
414
415 for( reactionIndex = 0; reactionIndex < numberOfReactions( ); ++reactionIndex ) {
416 if( m_reactionCrossSections[reactionIndex]->threshold( ) < URR_domainMax( ) ) {
417 reactions_in_URR_region.push_back( reactionIndex );
418 }
419 }
420
421 m_reactionsInURR_region.resize( reactions_in_URR_region.size( ) );
422 for( std::size_t i1 = 0; i1 < reactions_in_URR_region.size( ); ++i1 ) m_reactionsInURR_region[i1] = reactions_in_URR_region[i1];
423
424 if( a_settings._URR_mode( ) == Transporting::URR_mode::ACE_URR_probabilityTables ) {
425 auto URR_iter = a_setupInfo.m_ACE_URR_probabilityTablesFromGIDI.find( URR_label );
426 if( URR_iter != a_setupInfo.m_ACE_URR_probabilityTablesFromGIDI.end( ) ) {
427 auto ACE_URR_probabilityTablesIter = (*URR_iter).second->m_ACE_URR_probabilityTables.find( "total" );
428 if( ACE_URR_probabilityTablesIter != (*URR_iter).second->m_ACE_URR_probabilityTables.end( ) ) {
429 m_ACE_URR_probabilityTables = (*ACE_URR_probabilityTablesIter).second;
430 (*ACE_URR_probabilityTablesIter).second = nullptr; // Set to nullptr so destructor of m_ACE_URR_probabilityTables does not delete.
431 }
432 }
433 }
434 }
435
436 if( a_settings.addExpectedValueData( ) ) {
437 m_depositionEnergy.resize( totalCrossSection.length( ), 0.0 );
438 m_depositionMomentum.resize( totalCrossSection.length( ), 0.0 );
439 m_productionEnergy.resize( totalCrossSection.length( ), 0.0 );
440
441 m_gains.resize( a_particles.particles( ).size( ) );
442 std::size_t i2 = 0;
443 int projectileGainIndex = -1;
444 int photonGainIndex = -1;
445 for( std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( );
446 ++particle, ++i2 ) {
447 int particleIntid = a_setupInfo.m_particleIntids[particle->first];
448 int particleIndex = a_setupInfo.m_particleIndices[particle->first];
449
450 if( particleIntid == a_setupInfo.m_protare.projectileIntid( ) ) projectileGainIndex = (int) i2;
451 m_gains[i2] = new ContinuousEnergyGain( particleIntid, particleIndex, totalCrossSection.length( ) );
452 if( particle->first == PoPI::IDs::photon ) photonGainIndex = (int) i2;
453 }
454
455 std::vector< std::vector<double> > gains( a_particles.particles( ).size( ) );
456 for( std::size_t reactionIndex2 = 0; reactionIndex2 < a_reactions.size( ) + a_orphanProducts.size( ); ++reactionIndex2 ) {
457
458 std::size_t offset = 0;
459 std::vector<double> deposition_energy( m_energies.size( ), 0.0 );
460 std::vector<double> deposition_momentum( m_energies.size( ), 0.0 );
461 std::vector<double> production_energy( m_energies.size( ), 0.0 );
462 for( std::size_t i1 = 0; i1 < gains.size( ); ++i1 ) gains[i1] = std::vector<double>( m_energies.size( ), 0.0 );
463
464 HeatedReactionCrossSectionContinuousEnergy *MCGIDI_reaction_cross_section = nullptr;
465 GIDI::Reaction const *reaction = nullptr;
466 GIDI::Functions::Ys1d const *reactionCrossSection = nullptr; // If nullptr, reaction is an orphanProduct node.
467 GIDI::Functions::Function1dForm const *available_energy = nullptr;
468 GIDI::Functions::Function1dForm const *available_momentum = nullptr;
469 if( reactionIndex2 < a_reactions.size( ) ) {
470 reaction = a_reactions[reactionIndex2];
471 available_energy = reaction->availableEnergy( ).get<GIDI::Functions::Function1dForm>( 0 );
472 available_momentum = reaction->availableMomentum( ).get<GIDI::Functions::Function1dForm>( 0 );
473 MCGIDI_reaction_cross_section = m_reactionCrossSections[reactionIndex2];
474 offset = MCGIDI_reaction_cross_section->offset( ); }
475 else {
476 reaction = a_orphanProducts[reactionIndex2-a_reactions.size( )];
477 GIDI::Suite const &reactionCrossSectionSuite = reaction->crossSection( );
478 reactionCrossSection = reactionCrossSectionSuite.get<GIDI::Functions::Ys1d>( label );
479 offset = reactionCrossSection->start( );
480 }
481 if( a_settings.useSlowerContinuousEnergyConversion( ) ) { // Old way which is slow as it does one energy at a time.
482 for( std::size_t energy_index = offset; energy_index < m_energies.size( ); ++energy_index ) {
483 double energy = m_energies[energy_index];
484
485 if( reactionCrossSection == nullptr ) {
486 if( isPhotoAtomic ) { // Treat as Q = 0.0 since 2 photons will be emitted.
487 deposition_energy[energy_index] = energy; }
488 else {
489 deposition_energy[energy_index] = available_energy->evaluate( energy );
490
491 double Q = deposition_energy[energy_index] - energy; // Should use Q node to get this.
492 if( fabs( Q ) < 1e-12 * deposition_energy[energy_index] )
493 Q = 0.0; // Probably 0.0 due to rounding errors.
494 production_energy[energy_index] = Q;
495 }
496 deposition_momentum[energy_index] = available_momentum->evaluate( energy );
497 }
498
499 std::size_t i1 = 0;
500 for( std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( );
501 ++particle, ++i1 ) {
502 double product_energy, product_momentum, product_gain;
503
504 if( particle->first == PoPI::IDs::electron ) continue; // As of this coding, electrons are not complete in GNDS files.
505 // When they are, this statement can be removed.
506 if( ( reactionCrossSection != nullptr ) && ( particle->first != PoPI::IDs::photon ) ) continue;
507
508 if( reaction->isPairProduction( ) && ( particle->first == PoPI::IDs::photon ) ) {
509 product_energy = 2.0 * PoPI_electronMass_MeV_c2; // Assumes energy unit is MeV.
510 product_momentum = 0.0;
511 product_gain = 2.0; }
512 else {
513 reaction->continuousEnergyProductData( a_settings, particle->first, energy, product_energy, product_momentum,
514 product_gain, true );
515 }
516 if( (int) i1 == projectileGainIndex ) --product_gain;
517
518 deposition_energy[energy_index] -= product_energy;
519 deposition_momentum[energy_index] -= product_momentum;
520 gains[i1][energy_index] = product_gain;
521 }
522 } }
523 else { // New way which is hopefully faster.
524 if( reactionCrossSection == nullptr ) {
525 if( isPhotoAtomic ) { // Treat as Q = 0.0 since 2 photons will be emitted.
526 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) {
527 deposition_energy[energyIndex] = m_energies[energyIndex];
528 } }
529 else {
530 available_energy->mapToXsAndAdd( offset, *energiesPointer, deposition_energy, 1.0 );
531 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) {
532 double Q = deposition_energy[energyIndex] - m_energies[energyIndex];
533 if( fabs( Q ) < 1e-12 * deposition_energy[energyIndex] )
534 Q = 0.0; // Probably 0.0 due to rounding errors.
535 production_energy[energyIndex] = Q;
536 }
537 }
538 available_momentum->mapToXsAndAdd( offset, *energiesPointer, deposition_momentum, 1.0 );
539 }
540
541 std::size_t i1 = 0;
542 for( std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = a_particles.particles( ).begin( );
543 particle != a_particles.particles( ).end( ); ++particle, ++i1 ) {
544 if( particle->first == PoPI::IDs::electron ) continue; // As of this coding, electrons are not complete in GNDS files.
545 // When they are, this statement can be removed.
546 if( ( reactionCrossSection != nullptr ) && ( particle->first != PoPI::IDs::photon ) ) continue;
547
548 if( reaction->isPairProduction( ) && ( particle->first == PoPI::IDs::photon ) ) {
549 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) {
550 deposition_energy[energyIndex] -= 2.0 * PoPI_electronMass_MeV_c2; // Assumes energy unit is MeV.
551 gains[i1][energyIndex] = 2.0;
552 } }
553 else {
554 reaction->mapContinuousEnergyProductData( a_settings, particle->first, *energiesPointer, offset, deposition_energy,
555 deposition_momentum, gains[i1], true );
556 }
557
558 if( (int) i1 == projectileGainIndex ) {
559 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) --gains[i1][energyIndex];
560 }
561 }
562 }
563
564 if( a_particles.hasParticle( PoPI::IDs::photon ) ) {
565 if( a_setupInfo.m_initialStateIndices.find( reaction->label( ) ) != a_setupInfo.m_initialStateIndices.end( ) ) {
566 int intInitialStateIndex = a_setupInfo.m_initialStateIndices[reaction->label( )];
567 if( intInitialStateIndex >= 0 ) {
568 std::size_t initialStateIndex = static_cast<std::size_t>( intInitialStateIndex );
569 NuclideGammaBranchStateInfo *nuclideGammaBranchStateInfo = a_setupInfo.m_protare.nuclideGammaBranchStateInfos( )[initialStateIndex];
570 double multiplicity = nuclideGammaBranchStateInfo->multiplicity( );
571 double averageGammaEnergy = nuclideGammaBranchStateInfo->averageGammaEnergy( );
572
573 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) {
574 deposition_energy[energyIndex] -= averageGammaEnergy;
575 if( photonGainIndex >= 0 ) gains[static_cast<std::size_t>(photonGainIndex)][energyIndex] += multiplicity;
576 }
577 }
578 }
579 }
580
581 double crossSection = 0.0;
582 for( std::size_t energyIndex = offset; energyIndex < m_energies.size( ); ++energyIndex ) {
583 if( reactionCrossSection == nullptr ) {
584 crossSection = MCGIDI_reaction_cross_section->crossSection( energyIndex ); }
585 else {
586 crossSection = (*reactionCrossSection)[energyIndex-offset];
587 }
588
589 m_depositionEnergy[energyIndex] += crossSection * deposition_energy[energyIndex];
590 m_depositionMomentum[energyIndex] += crossSection * deposition_momentum[energyIndex];
591 m_productionEnergy[energyIndex] += crossSection * production_energy[energyIndex];
592 for( std::size_t i1 = 0; i1 < m_gains.size( ); ++i1 ) {
593 m_gains[i1]->adjustGain( energyIndex, crossSection * gains[i1][energyIndex] );
594 }
595 }
596 }
597 }
598}
G4ThreadLocal T * G4GeomSplitter< T >::offset
#define PoPI_electronMass_MeV_c2
Definition PoPI.hpp:31
virtual void mapToXsAndAdd(std::size_t a_offset, std::vector< double > const &a_Xs, std::vector< double > &a_results, double a_scaleFactor) const
Definition GIDI_form.cc:424
virtual double evaluate(double a_x1) const =0
void set(std::size_t a_index, double a_value)
Definition GIDI.hpp:1170
std::size_t start() const
Definition GIDI.hpp:1163
Grid const & grid() const
Definition GIDI.hpp:3368
bool has(std::string const &a_label) const
Definition GIDI.hpp:2616
T * get(std::size_t a_Index)
Definition GIDI.hpp:2642
LUPI_HOST_DEVICE Vector< double > const & energies() const
Definition MCGIDI.hpp:678
LUPI_HOST_DEVICE double threshold(std::size_t a_reactionIndex) const
Definition MCGIDI.hpp:651
LUPI_HOST_DEVICE Vector< MCGIDI_FLOAT > & totalCrossSection()
Definition MCGIDI.hpp:658
LUPI_HOST HeatedReactionCrossSectionContinuousEnergy const * reactionCrossSection(std::size_t a_index) const
Definition MCGIDI.hpp:640
LUPI_HOST_DEVICE std::size_t numberOfReactions() const
Definition MCGIDI.hpp:646
LUPI_HOST_DEVICE double crossSection(URR_protareInfos const &a_URR_protareInfos, int a_URR_index, std::size_t a_hashIndex, double a_energy, bool a_sampling=false) const
G4double energy(const ThreeVector &p, const G4double m)
LUPI_HOST ProbabilityBase2d * parseProbability2d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite, SetupInfo *a_setupInfo)
LUPI_HOST_DEVICE int binarySearchVector(double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
Definition MCGIDI.hpp:318
@ ptwXY_interpolationLinLin
Definition ptwXY.h:37
static std::string const photon
Definition PoPI.hpp:162
static std::string const electron
Definition PoPI.hpp:163

◆ ~HeatedCrossSectionContinuousEnergy()

LUPI_HOST_DEVICE MCGIDI::HeatedCrossSectionContinuousEnergy::~HeatedCrossSectionContinuousEnergy ( )

Definition at line 603 of file MCGIDI_heatedCrossSections.cc.

603 {
604
605 for( auto iter = m_reactionCrossSections.begin( ); iter < m_reactionCrossSections.end( ); ++iter ) delete *iter;
606 for( auto iter = m_gains.begin( ); iter < m_gains.end( ); ++iter ) delete *iter;
607 delete m_ACE_URR_probabilityTables;
608}

Member Function Documentation

◆ crossSection()

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

Definition at line 682 of file MCGIDI_heatedCrossSections.cc.

683 {
684
685 double energy_fraction;
686 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
687
688 if( a_URR_index >= 0 ) {
689 URR_protareInfo const &URR_protare_info = a_URR_protareInfos[static_cast<std::size_t>(a_URR_index)];
690
691 if( URR_protare_info.m_inURR ) {
692 double cross_section = 0.0;
693 for( std::size_t i1 = 0; i1 < m_reactionsInURR_region.size( ); ++i1 ) {
694 cross_section += reactionCrossSection2( m_reactionsInURR_region[i1], a_URR_protareInfos, a_URR_index, a_energy,
695 energy_index, energy_fraction, false );
696 }
697 return( cross_section );
698 }
699 }
700
701 return( energy_fraction * m_totalCrossSection[energy_index] + ( 1.0 - energy_fraction ) * m_totalCrossSection[energy_index+1] );
702}
LUPI_HOST_DEVICE double reactionCrossSection2(std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, int a_URR_index, double a_energy, std::size_t a_energyIndex, double a_energyFraction, bool a_sampling=false) const
LUPI_HOST_DEVICE std::size_t evaluationInfo(std::size_t a_hashIndex, double a_energy, double *a_energyFraction) const

Referenced by crossSectionAsGIDI_XYs1d(), and HeatedCrossSectionContinuousEnergy().

◆ crossSectionAsGIDI_XYs1d()

LUPI_HOST GIDI::Functions::XYs1d MCGIDI::HeatedCrossSectionContinuousEnergy::crossSectionAsGIDI_XYs1d ( ) const

Returns the total cross section as a GIDI::Functions::XYs1d instance.

Returns
A GIDI::Functions::XYs1d instance.

Definition at line 770 of file MCGIDI_heatedCrossSections.cc.

770 {
771
772 std::vector<double> energies = vectorToSTD_vector( m_energies );
773 std::vector<double> crossSection = vectorToSTD_vector( m_totalCrossSection );
774
775 return( GIDI::Functions::XYs1d( GIDI::Axes( ), ptwXY_interpolationLinLin, energies, crossSection, 0, m_temperature ) );
776}
LUPI_HOST std::vector< double > vectorToSTD_vector(Vector< double > a_input)

Referenced by reactionCrossSectionAsGIDI_XYs1d().

◆ depositionEnergy()

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

Returns the index of a sampled reaction for 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 projectile energy hash index.
a_energy[in] The energy of the projectile.

Definition at line 799 of file MCGIDI_heatedCrossSections.cc.

799 {
800
801 double energy_fraction;
802 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
803
804 return( energy_fraction * m_depositionEnergy[energy_index] + ( 1.0 - energy_fraction ) * m_depositionEnergy[energy_index+1] );
805}

◆ depositionMomentum()

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

Returns the index of a sampled reaction for 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 projectile energy hash index.
a_energy[in] The energy of the projectile.

Definition at line 815 of file MCGIDI_heatedCrossSections.cc.

815 {
816
817 double energy_fraction;
818 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
819
820 return( energy_fraction * m_depositionMomentum[energy_index] + ( 1.0 - energy_fraction ) * m_depositionMomentum[energy_index+1] );
821}

◆ energies()

LUPI_HOST_DEVICE Vector< double > const & MCGIDI::HeatedCrossSectionContinuousEnergy::energies ( ) const
inline

Returns a reference to m_styles.

Definition at line 678 of file MCGIDI.hpp.

Referenced by crossSectionAsGIDI_XYs1d(), and HeatedCrossSectionContinuousEnergy().

◆ evaluationInfo()

LUPI_HOST_DEVICE std::size_t MCGIDI::HeatedCrossSectionContinuousEnergy::evaluationInfo ( std::size_t a_hashIndex,
double a_energy,
double * a_energyFraction ) const

This function returns the index in a_energies where a_energy lies between the returned index and the next index. The returned index must lie between a_hashIndices[a_hashIndex] and a_hashIndices[a_hashIndex+1]. If a_energy is below the domain of a_energies, 0 is returned. If a_energy is above the domain of a_energies, the size of a_energies minus 2 is returned. The argument a_energyFraction the weight for the energy at the returned index with the next index getting weighting 1 minus a_energyFraction.

Parameters
a_hashIndex[in] Specifies projectile energy hash index.
a_energy[in] The energy whose index is requested.
a_energyFraction[in] This represents the weighting to apply to the two bounding energies.
Returns
The index bounding a_energy in the member m_energies.

Definition at line 626 of file MCGIDI_heatedCrossSections.cc.

626 {
627
628 return Sampling::evaluationForHashIndex( a_hashIndex, m_hashIndices, a_energy, m_energies, a_energyFraction );
629}
LUPI_HOST_DEVICE std::size_t evaluationForHashIndex(std::size_t a_hashIndex, Vector< std::size_t > const &a_hashIndices, double a_energy, Vector< double > const &a_energies, double *a_energyFraction)

Referenced by crossSection(), depositionEnergy(), depositionMomentum(), gain(), gainViaIntid(), productionEnergy(), and MCGIDI::HeatedCrossSectionsContinuousEnergy::sampleReaction().

◆ gain()

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

Returns the gain for the particle with index a_particleIndex for projectile energy a_energy.

Parameters
a_hashIndex[in] Specifies projectile energy hash index.
a_energy[in] The energy of the projectile.
a_particleIndex[in] The index of the particle whose gain is requested.

Definition at line 847 of file MCGIDI_heatedCrossSections.cc.

847 {
848
849 double energy_fraction;
850 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
851
852 for( std::size_t i1 = 0; i1 < m_gains.size( ); ++i1 ) {
853 if( a_particleIndex == m_gains[i1]->particleIndex( ) ) return( m_gains[i1]->gain( energy_index, energy_fraction ) );
854 }
855
856 return( 0.0 );
857}
LUPI_HOST_DEVICE double gain(std::size_t a_hashIndex, double a_energy, int a_particleIndex) const

Referenced by gain(), and gainViaIntid().

◆ gainViaIntid()

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

Returns the gain for the particle with intid a_particleIntid for projectile energy a_energy.

Parameters
a_hashIndex[in] Specifies projectile energy hash index.
a_energy[in] The energy of the projectile.
a_particleIntid[in] The intid of the particle whose gain is requested.

Definition at line 867 of file MCGIDI_heatedCrossSections.cc.

867 {
868
869 double energy_fraction;
870 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
871
872 for( std::size_t i1 = 0; i1 < m_gains.size( ); ++i1 ) {
873 if( a_particleIntid == m_gains[i1]->particleIntid( ) ) return( m_gains[i1]->gain( energy_index, energy_fraction ) );
874 }
875
876 return( 0.0 );
877}

◆ hasURR_probabilityTables()

LUPI_HOST_DEVICE bool MCGIDI::HeatedCrossSectionContinuousEnergy::hasURR_probabilityTables ( ) const

Returns true if this has a unresolved resonance region (URR) data and false otherwise.

Returns
true is if this has a URR data.

Definition at line 637 of file MCGIDI_heatedCrossSections.cc.

637 {
638
639 return( m_URR_mode != Transporting::URR_mode::none );
640}

Referenced by HeatedCrossSectionContinuousEnergy(), and reactionHasURR_probabilityTables().

◆ maximumEnergy()

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

Returns the maximum cross section domain.

Definition at line 645 of file MCGIDI.hpp.

◆ minimumEnergy()

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

Returns the minimum cross section domain.

Definition at line 644 of file MCGIDI.hpp.

◆ numberOfReactions()

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

Returns the number of reaction cross section.

Definition at line 646 of file MCGIDI.hpp.

Referenced by HeatedCrossSectionContinuousEnergy().

◆ print()

LUPI_HOST void MCGIDI::HeatedCrossSectionContinuousEnergy::print ( ProtareSingle const * a_protareSingle,
std::string const & a_indent,
std::string const & a_iFormat,
std::string const & a_energyFormat,
std::string const & a_dFormat ) const

Print to std::cout the content of this. This is mainly meant for debugging.

Parameters
a_protareSingle[in] The ProtareSingle instance this resides in.
a_indent[in] The buffer to read or write data to depending on a_mode.
a_iFormat[in] C printf format specifier for any interger that is printed (e.g., "%3d").
a_energyFormat[in] C printf format specifier for any interger that is printed (e.g., "%20.12e").
a_dFormat[in] C printf format specifier for any interger that is printed (e.g., "%14.7e").

Definition at line 979 of file MCGIDI_heatedCrossSections.cc.

980 {
981
982 char const *dFormat = a_dFormat.c_str( );
983 std::string indent2 = a_indent + " ";
984 std::string lineFormat = indent2 + a_energyFormat + " " + a_dFormat + " " + a_dFormat + " " + a_dFormat + " " + a_dFormat;
985
986 std::cout << std::endl;
987 std::cout << a_indent << "# Temperature = " << LUPI::Misc::doubleToString3( dFormat, m_temperature, true ) << std::endl;
988
989 std::cout << a_indent << "# Number of hash indices = " << m_hashIndices.size( ) << std::endl;
990 for( auto iter = m_hashIndices.begin( ); iter != m_hashIndices.end( ); ++iter )
991 std::cout << indent2 << LUPI::Misc::argumentsToString( a_iFormat.c_str( ), *iter ) << std::endl;
992 std::cout << std::endl;
993 std::size_t energySize = m_energies.size( );
994 std::cout << a_indent << "# Number of energies = " << m_energies.size( ) << std::endl;
995 std::cout << "# projectile total deposition deposition production" << std::endl;
996 std::cout << "# energy cross section energy momentum energy"<< std::endl;
997 std::cout << "#====================================================================================================="<< std::endl;
998 for( std::size_t index = 0; index != energySize; ++index ) {
999 std::cout << LUPI::Misc::argumentsToString( lineFormat.c_str( ), m_energies[index], m_totalCrossSection[index], m_depositionEnergy[index],
1000 m_depositionMomentum[index], m_productionEnergy[index] ) << std::endl;
1001 }
1002
1003 for( auto iter = m_gains.begin( ); iter != m_gains.end( ); ++iter ) (*iter)->print( a_protareSingle, a_indent, a_iFormat, a_energyFormat, a_dFormat );
1004
1005 std::cout << std::endl;
1006 std::cout << "# URR mode is ";
1007 if( m_URR_mode == Transporting::URR_mode::none ) {
1008 std::cout << "none" << std::endl; }
1009 else if( m_URR_mode == Transporting::URR_mode::pdfs ) {
1010 std::cout << "pdfs" << std::endl; }
1011 else if( m_URR_mode == Transporting::URR_mode::ACE_URR_probabilityTables ) {
1012 std::cout << "ACE style protability tables" << std::endl; }
1013 else {
1014 std::cout << "Oops, need to code into print method." << std::endl;
1015 }
1016
1017 std::cout << a_indent << "# Index of reactions in URR region:";
1018 for( auto reactionIter = m_reactionsInURR_region.begin( ); reactionIter != m_reactionsInURR_region.end( ); ++reactionIter ) {
1019 std::cout << indent2 << LUPI::Misc::argumentsToString( a_iFormat.c_str( ), *reactionIter ) << std::endl;
1020 }
1021 std::cout << std::endl;
1022
1023 std::size_t reactionIndex = 0;
1024 std::cout << std::endl;
1025 std::cout << a_indent << "# Number of reactions = " << m_reactionCrossSections.size( ) << std::endl;
1026 for( auto iter = m_reactionCrossSections.begin( ); iter != m_reactionCrossSections.end( ); ++iter, ++reactionIndex ) {
1027 Reaction const *reaction = a_protareSingle->reaction( reactionIndex );
1028
1029 std::cout << a_indent << "# Reaction number " << reactionIndex << std::endl;
1030 std::cout << a_indent << "# Reaction label: " << reaction->label( ).c_str( ) << std::endl;
1031 (*iter)->print( a_protareSingle, a_indent, a_iFormat, a_energyFormat, a_dFormat );
1032 }
1033}
std::string argumentsToString(char const *a_format,...)
Definition LUPI_misc.cc:305
std::string doubleToString3(char const *a_format, double a_value, bool a_reduceBits=false)
Definition LUPI_misc.cc:327

◆ productionEnergy()

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

Returns the index of a sampled reaction for 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 projectile energy hash index.
a_energy[in] The energy of the projectile.

Definition at line 831 of file MCGIDI_heatedCrossSections.cc.

831 {
832
833 double energy_fraction;
834 std::size_t energy_index = evaluationInfo( a_hashIndex, a_energy, &energy_fraction );
835
836 return( energy_fraction * m_productionEnergy[energy_index] + ( 1.0 - energy_fraction ) * m_productionEnergy[energy_index+1] );
837}

◆ reactionCrossSection() [1/3]

LUPI_HOST HeatedReactionCrossSectionContinuousEnergy const * MCGIDI::HeatedCrossSectionContinuousEnergy::reactionCrossSection ( std::size_t a_index) const
inline

Returns the reaction cross section at index a_index.

Definition at line 640 of file MCGIDI.hpp.

Referenced by HeatedCrossSectionContinuousEnergy(), URR_domainMax(), and URR_domainMin().

◆ reactionCrossSection() [2/3]

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::reactionCrossSection ( std::size_t a_reactionIndex,
URR_protareInfos const & a_URR_protareInfos,
int a_URR_index,
double a_energy ) const

Definition at line 742 of file MCGIDI_heatedCrossSections.cc.

743 {
744
745 int intEnergyIndex = binarySearchVector( a_energy_in, m_energies );
746 std::size_t energyIndex = static_cast<std::size_t>( intEnergyIndex );
747 double energyFraction;
748
749 if( intEnergyIndex < 0 ) {
750 if( intEnergyIndex == -1 ) {
751 energyIndex = m_energies.size( ) - 2;
752 energyFraction = 0.0; }
753 else {
754 energyIndex = 0;
755 energyFraction = 1.0;
756 } }
757 else {
758 energyFraction = ( m_energies[energyIndex+1] - a_energy_in ) / ( m_energies[energyIndex+1] - m_energies[energyIndex] );
759 }
760
761 return( reactionCrossSection2( a_reactionIndex, a_URR_protareInfos, a_URR_index, a_energy_in, energyIndex, energyFraction, false ) );
762}

◆ reactionCrossSection() [3/3]

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

◆ reactionCrossSection2()

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::reactionCrossSection2 ( std::size_t a_reactionIndex,
URR_protareInfos const & a_URR_protareInfos,
int a_URR_index,
double a_energy,
std::size_t a_energyIndex,
double a_energyFraction,
bool a_sampling = false ) const

Definition at line 717 of file MCGIDI_heatedCrossSections.cc.

718 {
719
720 HeatedReactionCrossSectionContinuousEnergy const &reaction = *m_reactionCrossSections[a_reactionIndex];
721 double URR_cross_section_factor = 1.0;
722
723 if( a_URR_index >= 0 ) {
724 URR_protareInfo const &URR_protare_info = a_URR_protareInfos[static_cast<std::size_t>(a_URR_index)];
725 if( URR_protare_info.m_inURR ) {
726 if( m_URR_mode == Transporting::URR_mode::pdfs ) {
727 if( reaction.URR_probabilityTables( ) != nullptr )
728 URR_cross_section_factor = reaction.URR_probabilityTables( )->sample( a_energy, URR_protare_info.m_rng_Value, []() -> double { return 0.0; } ); }
730 if( reaction._ACE_URR_probabilityTables( ) != nullptr )
731 URR_cross_section_factor = reaction._ACE_URR_probabilityTables( )->sample( a_energy, URR_protare_info.m_rng_Value );
732 }
733 }
734 }
735
736 return( URR_cross_section_factor * ( a_energyFraction * reaction.crossSection( a_energyIndex ) + ( 1.0 - a_energyFraction ) * reaction.crossSection( a_energyIndex + 1 ) ) );
737}

Referenced by crossSection(), reactionCrossSection(), and MCGIDI::HeatedCrossSectionsContinuousEnergy::sampleReaction().

◆ reactionCrossSectionAsGIDI_XYs1d()

LUPI_HOST GIDI::Functions::XYs1d MCGIDI::HeatedCrossSectionContinuousEnergy::reactionCrossSectionAsGIDI_XYs1d ( std::size_t a_reactionIndex) const

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

Parameters
a_reactionIndex[in] Specifies the indexs of the reaction whose cross section is requested.
Returns
A GIDI::Functions::XYs1d instance.

Definition at line 786 of file MCGIDI_heatedCrossSections.cc.

786 {
787
788 return( m_reactionCrossSections[a_reactionIndex]->crossSectionAsGIDI_XYs1d( m_temperature, m_energies ) );
789}
LUPI_HOST GIDI::Functions::XYs1d crossSectionAsGIDI_XYs1d() const

◆ reactionHasURR_probabilityTables()

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

Definition at line 656 of file MCGIDI.hpp.

656{ return( m_reactionCrossSections[a_index]->hasURR_probabilityTables( ) ); }

◆ serialize()

LUPI_HOST_DEVICE void MCGIDI::HeatedCrossSectionContinuousEnergy::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 911 of file MCGIDI_heatedCrossSections.cc.

911 {
912
913 DATA_MEMBER_DOUBLE( m_temperature, a_buffer, a_mode );
914 DATA_MEMBER_VECTOR_SIZE_T( m_hashIndices, a_buffer, a_mode );
915 DATA_MEMBER_VECTOR_DOUBLE( m_energies, a_buffer, a_mode );
916 DATA_MEMBER_VECTOR_FLOAT_OR_DOUBLE( m_totalCrossSection, a_buffer, a_mode );
917 DATA_MEMBER_VECTOR_FLOAT_OR_DOUBLE( m_depositionEnergy, a_buffer, a_mode );
918 DATA_MEMBER_VECTOR_FLOAT_OR_DOUBLE( m_depositionMomentum, a_buffer, a_mode );
919 DATA_MEMBER_VECTOR_FLOAT_OR_DOUBLE( m_productionEnergy, a_buffer, a_mode );
920 m_URR_mode = serializeURR_mode( m_URR_mode, a_buffer, a_mode );
921 DATA_MEMBER_VECTOR_SIZE_T( m_reactionsInURR_region, a_buffer, a_mode );
922 m_ACE_URR_probabilityTables = serializeACE_URR_probabilityTables( m_ACE_URR_probabilityTables, a_buffer, a_mode );
923
924 std::size_t vectorSize = m_reactionCrossSections.size( );
925 int vectorSizeInt = (int) vectorSize;
926 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
927 vectorSize = (std::size_t) vectorSizeInt;
928
929 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) m_reactionCrossSections.resize( vectorSize, &a_buffer.m_placement );
930 if( a_mode == LUPI::DataBuffer::Mode::Memory ) a_buffer.m_placement += m_reactionCrossSections.internalSize();
931 for( std::size_t memberIndex = 0; memberIndex < vectorSize; ++memberIndex ) {
932 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
933 if( a_buffer.m_placement != nullptr ) {
934 m_reactionCrossSections[memberIndex] = new(a_buffer.m_placement) HeatedReactionCrossSectionContinuousEnergy;
935 a_buffer.incrementPlacement( sizeof( HeatedReactionCrossSectionContinuousEnergy ) ); }
936 else {
937 m_reactionCrossSections[memberIndex] = new HeatedReactionCrossSectionContinuousEnergy;
938 }
939 }
940 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
941 a_buffer.incrementPlacement( sizeof( HeatedReactionCrossSectionContinuousEnergy ) );
942 }
943 m_reactionCrossSections[memberIndex]->serialize( a_buffer, a_mode );
944 }
945
946 vectorSize = m_gains.size( );
947 vectorSizeInt = (int) vectorSize;
948 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
949 vectorSize = (std::size_t) vectorSizeInt;
950 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) m_gains.resize( vectorSize, &a_buffer.m_placement );
951 if( a_mode == LUPI::DataBuffer::Mode::Memory ) a_buffer.m_placement += m_gains.internalSize();
952 for( std::size_t memberIndex = 0; memberIndex < vectorSize; ++memberIndex ) {
953 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
954 if( a_buffer.m_placement != nullptr ) {
955 m_gains[memberIndex] = new(a_buffer.m_placement) ContinuousEnergyGain;
956 a_buffer.incrementPlacement( sizeof( ContinuousEnergyGain ) ); }
957 else {
958 m_gains[memberIndex] = new ContinuousEnergyGain;
959 }
960 }
961
962 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
963 a_buffer.incrementPlacement( sizeof( ContinuousEnergyGain ) );
964 }
965 m_gains[memberIndex]->serialize( a_buffer, a_mode );
966 }
967}
#define DATA_MEMBER_VECTOR_DOUBLE(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)
#define DATA_MEMBER_VECTOR_FLOAT_OR_DOUBLE
Definition MCGIDI.hpp:19
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
LUPI_HOST_DEVICE ACE_URR_probabilityTables * serializeACE_URR_probabilityTables(ACE_URR_probabilityTables *a_ACE_URR_probabilityTables, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Transporting::URR_mode serializeURR_mode(Transporting::URR_mode a_URR_mode, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

◆ setUserParticleIndex()

LUPI_HOST void MCGIDI::HeatedCrossSectionContinuousEnergy::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 886 of file MCGIDI_heatedCrossSections.cc.

886 {
887
888 for( auto iter = m_gains.begin( ); iter != m_gains.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
889}

◆ setUserParticleIndexViaIntid()

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

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

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

Definition at line 898 of file MCGIDI_heatedCrossSections.cc.

898 {
899
900 for( auto iter = m_gains.begin( ); iter != m_gains.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
901}

◆ temperature()

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::temperature ( ) const
inline

Returns the value of the m_temperature member.

Definition at line 643 of file MCGIDI.hpp.

Referenced by HeatedCrossSectionContinuousEnergy().

◆ threshold()

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::threshold ( std::size_t a_reactionIndex) const
inline

Returns the threshold for the reaction with index a_reactionIndex.

Definition at line 651 of file MCGIDI.hpp.

Referenced by HeatedCrossSectionContinuousEnergy(), and threshold().

◆ thresholdOffset()

LUPI_HOST_DEVICE std::size_t MCGIDI::HeatedCrossSectionContinuousEnergy::thresholdOffset ( std::size_t a_reactionIndex) const
inline

Returns the offset for the cross section for the reaction with index a_reactionIndex.

Definition at line 649 of file MCGIDI.hpp.

◆ totalCrossSection()

LUPI_HOST_DEVICE Vector< MCGIDI_FLOAT > & MCGIDI::HeatedCrossSectionContinuousEnergy::totalCrossSection ( )
inline

Returns a reference to member m_totalCrossSection.

Definition at line 658 of file MCGIDI.hpp.

Referenced by HeatedCrossSectionContinuousEnergy().

◆ URR_domainMax()

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::URR_domainMax ( ) const

Returns the maximum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.

Returns
true is if this has a URR data.

Definition at line 667 of file MCGIDI_heatedCrossSections.cc.

667 {
668
669 if( m_ACE_URR_probabilityTables != nullptr ) return( m_ACE_URR_probabilityTables->domainMax( ) );
670
671 for( std::size_t i1 = 0; i1 < m_reactionCrossSections.size( ); ++i1 ) {
672 HeatedReactionCrossSectionContinuousEnergy *reactionCrossSection = m_reactionCrossSections[i1];
673
674 if( reactionCrossSection->hasURR_probabilityTables( ) ) return( reactionCrossSection->URR_domainMax( ) );
675 }
676
677 return( -1.0 );
678}

Referenced by HeatedCrossSectionContinuousEnergy().

◆ URR_domainMin()

LUPI_HOST_DEVICE double MCGIDI::HeatedCrossSectionContinuousEnergy::URR_domainMin ( ) const

Returns the minimum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.

Returns
true is if this has a URR data.

Definition at line 648 of file MCGIDI_heatedCrossSections.cc.

648 {
649
650 if( m_ACE_URR_probabilityTables != nullptr ) return( m_ACE_URR_probabilityTables->domainMin( ) );
651
652 for( std::size_t i1 = 0; i1 < m_reactionCrossSections.size( ); ++i1 ) {
653 HeatedReactionCrossSectionContinuousEnergy *reactionCrossSection = m_reactionCrossSections[i1];
654
655 if( reactionCrossSection->hasURR_probabilityTables( ) ) return( reactionCrossSection->URR_domainMin( ) );
656 }
657
658 return( -1.0 );
659}

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