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

Simple C++ string class, useful as replacement for std::string if this cannot be used, or just for fun. More...

Namespaces

namespace  Transporting
namespace  Distributions
namespace  Functions
namespace  Probabilities
namespace  Sampling

Classes

class  ACE_URR_probabilityTablesFromGIDI
class  SetupInfo
class  MultiGroupHash
class  URR_protareInfo
class  URR_protareInfos
class  ACE_URR_probabilityTable
class  ACE_URR_probabilityTables
class  HeatedReactionCrossSectionContinuousEnergy
class  ContinuousEnergyGain
class  HeatedCrossSectionContinuousEnergy
class  HeatedCrossSectionsContinuousEnergy
class  MultiGroupGain
class  HeatedReactionCrossSectionMultiGroup
class  HeatedCrossSectionMultiGroup
class  HeatedCrossSectionsMultiGroup
class  NuclideGammaBranchInfo
class  NuclideGammaBranchStateInfo
class  GRIN_levelsAndProbabilities
class  GRIN_inelasticForEnergy
class  GRIN_inelastic
class  GRIN_captureToCompound
class  GRIN_captureLevelProbability
class  GRIN_capture
class  Product
class  DelayedNeutron
class  OutputChannel
class  Reaction
class  Protare
class  ProtareSingle
class  ProtareComposite
class  ProtareTNSL
class  DomainHash
class  String
class  Vector

Enumerations

enum class  ProtareType { single , composite , TNSL }
enum class  TwoBodyOrder { notApplicable , firstParticle , secondParticle }
enum class  ChannelType { none , twoBody , uncorrelatedBodies }
enum class  Interpolation {
  LINLIN , LINLOG , LOGLIN , LOGLOG ,
  FLAT , OTHER
}
enum class  Function1dType {
  none , constant , XYs , polyomial ,
  gridded , regions , branching , TerrellFissionNeutronMultiplicityModel
}
enum class  Function2dType { none , XYs }
enum class  ProbabilityBase1dType { none , xs_pdf_cdf }
enum class  ProbabilityBase2dType {
  none , XYs , regions , isotropic ,
  discreteGamma , primaryGamma , recoil , NBodyPhaseSpace ,
  evaporation , generalEvaporation , simpleMaxwellianFission , Watt ,
  weightedFunctionals
}
enum class  ProbabilityBase3dType { none , XYs }

Functions

LUPI_HOST_DEVICE double particleKineticEnergy (double a_mass_unitOfEnergy, double a_particleBeta)
LUPI_HOST_DEVICE double particleKineticEnergyFromBeta2 (double a_mass_unitOfEnergy, double a_particleBeta2)
LUPI_HOST_DEVICE double boostSpeed (double a_massProjectile, double a_kineticEnergyProjectile, double a_massTarget)
LUPI_HOST_DEVICE int muCOM_From_muLab (double a_muLab, double a_boostBeta, double a_productBeta, double &a_muPlus, double &a_JacobianPlus, double &a_muMinus, double &a_JacobianMinus)
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_DEVICE int binarySearchVector (double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
LUPI_HOST_DEVICE int binarySearchVectorBounded (double a_x, Vector< double > const &a_Xs, std::size_t a_lower, std::size_t a_upper, bool a_boundIndex)
LUPI_HOST ProtareprotareFromGIDIProtare (LUPI::StatusMessageReporting &a_smr, GIDI::Protare 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 Vector< double > GIDI_VectorDoublesToMCGIDI_VectorDoubles (GIDI::Vector a_vector)
LUPI_HOST void addVectorItemsToSet (Vector< int > const &a_from, std::set< int > &a_to)
LUPI_HOST_DEVICE int distributionTypeToInt (Distributions::Type a_type)
LUPI_HOST_DEVICE Distributions::Type intToDistributionType (int a_type)
LUPI_HOST_DEVICE void serializeProducts (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< Product * > &a_products)
LUPI_HOST_DEVICE void serializeDelayedNeutrons (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< DelayedNeutron * > &a_delayedNeutrons)
LUPI_HOST_DEVICE void serializeQs (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< Functions::Function1d_d1 * > &a_Qs)
LUPI_HOST_DEVICE void serializeFissionResiduals (GIDI::Construction::FissionResiduals &a_fissionResiduals, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST void convertACE_URR_probabilityTablesFromGIDI (GIDI::ProtareSingle const &a_protare, Transporting::MC &a_settings, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE Transporting::URR_mode serializeURR_mode (Transporting::URR_mode a_URR_mode, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE ACE_URR_probabilityTablesserializeACE_URR_probabilityTables (ACE_URR_probabilityTables *a_ACE_URR_probabilityTables, LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Distributions::DistributionserializeDistribution (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Distributions::Distribution *a_distribution)
LUPI_HOST std::vector< double > vectorToSTD_vector (Vector< double > a_input)
LUPI_HOST std::vector< double > vectorToSTD_vector (Vector< float > a_input)
LUPI_HOST_DEVICE Interpolation GIDI2MCGIDI_interpolation (ptwXY_interpolation a_interpolation)
LUPI_HOST_DEVICE Function1dType Function1dClass (Functions::Function1d *funct)
LUPI_HOST_DEVICE Functions::Function1dserializeFunction1d (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d *a_function1d)
LUPI_HOST_DEVICE Functions::Function1d_d1serializeFunction1d_d1 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d)
LUPI_HOST_DEVICE Functions::Function1d_d2serializeFunction1d_d2 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d2 *a_function1d)
LUPI_HOST_DEVICE Function2dType Function2dClass (Functions::Function2d *funct)
LUPI_HOST_DEVICE Functions::Function2dserializeFunction2d (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function2d *a_function2d)
LUPI_HOST_DEVICE ProbabilityBase1dType ProbabilityBase1dClass (Probabilities::ProbabilityBase1d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase1dserializeProbability1d (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase1d *a_probability1d)
LUPI_HOST_DEVICE ProbabilityBase2dType ProbabilityBase2dClass (Probabilities::ProbabilityBase2d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2dserializeProbability2d (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d *a_probability2d)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1serializeProbability2d_d1 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d1 *a_probability2d)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d2serializeProbability2d_d2 (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d2 *a_probability2d)
LUPI_HOST_DEVICE ProbabilityBase3dType ProbabilityBase3dClass (Probabilities::ProbabilityBase3d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3dserializeProbability3d (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase3d *a_probability3d)
template<typename RNG>
LUPI_HOST_DEVICE void upScatterModelABoostParticle (Sampling::Input &a_input, RNG &&a_rng, Sampling::Product &a_product)
LUPI_HOST_DEVICE bool operator< (const String &, const String &)
LUPI_HOST_DEVICE int MCGIDI_strcmp (const char *p1, const char *p2)
LUPI_HOST_DEVICE size_t MCGIDI_strlen (const char *str)
LUPI_HOST_DEVICE void MCGIDI_memmove (char *dest, const char *src, size_t n)
LUPI_HOST_DEVICE int MCGIDI_strncmp (const char *s1, const char *s2, size_t n)
LUPI_HOST_DEVICE String operator+ (const String &lhs, const String &rhs)

Detailed Description

Simple C++ string class, useful as replacement for std::string if this cannot be used, or just for fun.

Enumeration Type Documentation

◆ ChannelType

enum class MCGIDI::ChannelType
strong
Enumerator
none 
twoBody 
uncorrelatedBodies 

Definition at line 381 of file MCGIDI.hpp.

◆ Function1dType

enum class MCGIDI::Function1dType
strong
Enumerator
none 
constant 
XYs 
polyomial 
gridded 
regions 
branching 
TerrellFissionNeutronMultiplicityModel 

Definition at line 20 of file MCGIDI_functions.hpp.

◆ Function2dType

enum class MCGIDI::Function2dType
strong
Enumerator
none 
XYs 

Definition at line 21 of file MCGIDI_functions.hpp.

21{ none, XYs };

◆ Interpolation

enum class MCGIDI::Interpolation
strong
Enumerator
LINLIN 
LINLOG 
LOGLIN 
LOGLOG 
FLAT 
OTHER 

Definition at line 19 of file MCGIDI_functions.hpp.

◆ ProbabilityBase1dType

enum class MCGIDI::ProbabilityBase1dType
strong
Enumerator
none 
xs_pdf_cdf 

Definition at line 22 of file MCGIDI_functions.hpp.

◆ ProbabilityBase2dType

◆ ProbabilityBase3dType

enum class MCGIDI::ProbabilityBase3dType
strong
Enumerator
none 
XYs 

Definition at line 26 of file MCGIDI_functions.hpp.

26{ none, XYs };

◆ ProtareType

enum class MCGIDI::ProtareType
strong
Enumerator
single 
composite 
TNSL 

Definition at line 79 of file MCGIDI.hpp.

◆ TwoBodyOrder

enum class MCGIDI::TwoBodyOrder
strong
Enumerator
notApplicable 
firstParticle 
secondParticle 

Definition at line 229 of file MCGIDI.hpp.

Function Documentation

◆ addVectorItemsToSet()

LUPI_HOST void MCGIDI::addVectorItemsToSet ( Vector< int > const & a_from,
std::set< int > & a_to )

Adds the items in Vector a_from to the set a_to.

Parameters
a_to[in] The list of ints to add to the set.
a_from[in] The set to add the ints to.

Definition at line 113 of file MCGIDI_misc.cc.

113 {
114
115 for( Vector<int>::const_iterator iter = a_from.begin( ); iter != a_from.end( ); ++iter ) a_to.insert( *iter );
116}
LUPI_HOST_DEVICE iterator begin()
LUPI_HOST_DEVICE iterator end()

◆ binarySearchVector()

LUPI_HOST_DEVICE int MCGIDI::binarySearchVector ( double a_x,
Vector< double > const & a_Xs,
bool a_boundIndex = false )
inline

Definition at line 318 of file MCGIDI.hpp.

318 {
319
320 std::size_t lower = 0, middle, upper = a_Xs.size( ) - 1;
321
322 if( a_Xs.size( ) == 0 ) {
323 return( -3 ); }
324 else if( a_x < a_Xs[0] ) {
325 if( a_boundIndex ) return( 0 );
326 return( -2 ); }
327 else if( a_x > a_Xs.back( ) ) {
328 if( a_boundIndex ) return( static_cast<int>( upper ) );
329 return( -1 );
330 }
331
332 while( 1 ) {
333 middle = ( lower + upper ) >> 1;
334 if( middle == lower ) break;
335 if( a_x < a_Xs[middle] ) {
336 upper = middle; }
337 else {
338 lower = middle;
339 }
340 }
341 return( static_cast<int>( lower ) );
342}
LUPI_HOST_DEVICE std::size_t size() const
LUPI_HOST_DEVICE T & back()

Referenced by MCGIDI::Distributions::IncoherentBoundToFreePhotoAtomicScattering::angleBiasing(), MCGIDI::Distributions::CoherentPhotoAtomicScattering::evaluate(), MCGIDI::Functions::Gridded1d::evaluate(), MCGIDI::Functions::Regions1d::evaluate(), MCGIDI::Functions::XYs1d::evaluate(), MCGIDI::Functions::XYs2d::evaluate(), MCGIDI::Probabilities::Regions2d::evaluate(), MCGIDI::Probabilities::Xs_pdf_cdf1d::evaluate(), MCGIDI::Probabilities::XYs2d::evaluate(), MCGIDI::Probabilities::XYs3d::evaluate(), MCGIDI::Distributions::CoherentPhotoAtomicScattering::evaluateFormFactor(), MCGIDI::Distributions::IncoherentBoundToFreePhotoAtomicScattering::evaluateOccupationNumber(), MCGIDI::Distributions::IncoherentPhotoAtomicScattering::evaluateScatteringFactor(), MCGIDI::HeatedCrossSectionContinuousEnergy::HeatedCrossSectionContinuousEnergy(), MCGIDI::MultiGroupHash::index(), MCGIDI::HeatedCrossSectionContinuousEnergy::reactionCrossSection(), MCGIDI::HeatedCrossSectionsMultiGroup::reactionCrossSection(), MCGIDI::ACE_URR_probabilityTable::sample(), MCGIDI::ACE_URR_probabilityTables::sample(), MCGIDI::Distributions::CoherentElasticTNSL::sample(), MCGIDI::Distributions::CoherentPhotoAtomicScattering::sample(), MCGIDI::Distributions::IncoherentBoundToFreePhotoAtomicScattering::sample(), MCGIDI::Probabilities::Regions2d::sample(), MCGIDI::Probabilities::Xs_pdf_cdf1d::sample(), MCGIDI::Probabilities::XYs2d::sample(), MCGIDI::Probabilities::XYs3d::sample(), MCGIDI::Probabilities::XYs2d::sample2dOf3d(), and MCGIDI::ProtareSingle::sampleTargetBetaForUpscatterModelA().

◆ binarySearchVectorBounded()

LUPI_HOST_DEVICE int MCGIDI::binarySearchVectorBounded ( double a_x,
Vector< double > const & a_Xs,
std::size_t a_lower,
std::size_t a_upper,
bool a_boundIndex )
inline

Definition at line 347 of file MCGIDI.hpp.

348 {
349
350 std::size_t middle;
351
352 if( a_Xs.size( ) == 0 ) {
353 return( -3 ); }
354 else if( a_x < a_Xs[a_lower] ) {
355 if( a_boundIndex ) return( static_cast<int>( a_lower ) );
356 return( -2 ); }
357 else if( a_x > a_Xs[a_upper] ) {
358 if( a_boundIndex ) return( static_cast<int>( a_upper ) );
359 return( -1 );
360 }
361
362 while( 1 ) {
363 middle = ( a_lower + a_upper ) >> 1;
364 if( middle == a_lower ) break;
365 if( a_x < a_Xs[middle] ) {
366 a_upper = middle; }
367 else {
368 a_lower = middle;
369 }
370 }
371 return( static_cast<int>( a_lower ) );
372}

Referenced by MCGIDI::Sampling::evaluationForHashIndex().

◆ boostSpeed()

LUPI_HOST_DEVICE double MCGIDI::boostSpeed ( double a_massProjectile,
double a_kineticEnergyProjectile,
double a_massTarget )

This function returns the boost speed required to boost to the center-of-mass for a projectile hitting a target.

Parameters
a_massProjectile[in] The mass of the projectile in energy units.
a_kineticEnergyProjectile[in] The kinetic energy of the projectile.
a_massTarget[in] The mass of the target in energy units.
Returns
The relativistic kinetic energy of the particle.

Definition at line 161 of file MCGIDI_misc.cc.

161 {
162
163 double betaProjectile = MCGIDI_particleBeta( a_massProjectile, a_kineticEnergyProjectile );
164
165 return( betaProjectile / ( 1.0 + a_massTarget / ( a_massProjectile + a_kineticEnergyProjectile ) ) );
166}
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71

◆ convertACE_URR_probabilityTablesFromGIDI()

LUPI_HOST void MCGIDI::convertACE_URR_probabilityTablesFromGIDI ( GIDI::ProtareSingle const & a_protare,
Transporting::MC & a_settings,
SetupInfo & a_setupInfo )

This method serializes a ACE_URR_probabilityTables instance pointed to by a_ACE_URR_probabilityTables.

Parameters
a_protare[in] The GIDI::Protare whose data is to be used to construct this.
a_settings[in] Used to pass user options to the this to instruct it which data are desired.
a_setupInfo[in] Used internally when constructing a Protare to pass information to other constructors.
Returns
A pointer to the converted ACE_URR_probabilityTables instance.

Definition at line 282 of file MCGIDI_URR.cc.

282 {
283
286
287 int64_t numberConverted;
288 char *endCharacter;
289
290 for( auto iter = a_protare.ACE_URR_probabilityTables( ).begin( ); iter != a_protare.ACE_URR_probabilityTables( ).end( ); ++iter ) {
291 bool needToInitialize( true );
292 std::map<std::size_t, std::string> columnNames;
293 ACE_URR_probabilityTablesFromGIDI *ACE_URR_probabilityTablesFromGIDI1 = new ACE_URR_probabilityTablesFromGIDI( );
295 GIDI::ACE_URR::ProbabilityTable::Forms &incidentEnergies = form->forms( );
296
297 for( auto incidentEnergyIter = incidentEnergies.begin( ); incidentEnergyIter != incidentEnergies.end( ); ++incidentEnergyIter ) {
298 GIDI::ACE_URR::IncidentEnergy *incidentEnergy = *incidentEnergyIter;
299 GIDI::Table::Table const &table = incidentEnergy->table( );
300 std::size_t numberOfRows = static_cast<std::size_t>( table.rows( ) );
301 std::size_t numberOfColumns = static_cast<std::size_t>( table.columns( ) );
302
303 std::size_t columnIndex = 0;
304 for( auto columnHeaderIter = table.columnHeaders( ).begin( ); columnHeaderIter != table.columnHeaders( ).end( ); ++columnHeaderIter ) {
305 GIDI::Table::Column const *columnHeader = dynamic_cast<GIDI::Table::Column *>( *columnHeaderIter );
306 if( needToInitialize && columnIndex > 0 ) {
307 ACE_URR_probabilityTablesFromGIDI1->m_ACE_URR_probabilityTables[columnHeader->name( )] =
308 new ACE_URR_probabilityTables( incidentEnergies.size( ) );
309 columnNames[columnIndex] = columnHeader->name( );
310 }
311 ++columnIndex;
312 }
313 needToInitialize = false;
314
315 GIDI::Table::Data const &data = table.data( );
316
317 std::string const &body = data.body( );
318 char const *text = body.c_str( );
319 double *dValues = nfu_stringToListOfDoubles( NULL, text, data.sep( ).c_str( )[0], &numberConverted, &endCharacter, 0 );
320 if( dValues == nullptr ) throw GIDI::Exception( "convertACE_URR_probabilityTablesFromGIDI: nfu_stringToListOfDoubles failed." );
321
322 std::vector<std::vector<double> > columns( numberOfColumns );
323 for( columnIndex = 0; columnIndex < numberOfColumns; ++columnIndex ) {
324 columns[columnIndex].reserve( numberOfRows );
325 for( std::size_t rowIndex = 0; rowIndex < numberOfRows; ++rowIndex ) columns[columnIndex].push_back( dValues[rowIndex*numberOfColumns+columnIndex] );
326 }
327 free( dValues );
328
329 for( columnIndex = 1; columnIndex < numberOfColumns; ++columnIndex ) {
330 ACE_URR_probabilityTablesFromGIDI1->m_ACE_URR_probabilityTables[columnNames[columnIndex]]->push_back(
331 new ACE_URR_probabilityTable( incidentEnergy->value( ), columns[0], columns[columnIndex] ) );
332 }
333 }
334 a_setupInfo.m_ACE_URR_probabilityTablesFromGIDI[form->label()] = ACE_URR_probabilityTablesFromGIDI1;
335 }
336 }
337}
Table::Table const & table() const
Definition GIDI.hpp:4212
std::vector< IncidentEnergy * > Forms
Definition GIDI.hpp:4225
std::string const & label() const
Definition GIDI.hpp:658
std::string const & name() const
Definition GIDI.hpp:2795
std::string const & body() const
Definition GIDI.hpp:2821
std::string const & sep() const
Definition GIDI.hpp:2820
std::map< std::string, ACE_URR_probabilityTables * > m_ACE_URR_probabilityTables
Definition MCGIDI.hpp:243
std::map< std::string, ACE_URR_probabilityTablesFromGIDI * > m_ACE_URR_probabilityTablesFromGIDI
Definition MCGIDI.hpp:279
LUPI_HOST URR_mode _URR_mode() const
Definition MCGIDI.hpp:194
LUPI_HOST LookupMode::Data1d crossSectionLookupMode() const
Definition MCGIDI.hpp:158
void free(voidpf ptr)
double * nfu_stringToListOfDoubles(statusMessageReporting *smr, char const *str, char sep, int64_t *numberConverted, char **endCharacter, int useSystem_strtod)

◆ distributionTypeToInt()

LUPI_HOST_DEVICE int MCGIDI::distributionTypeToInt ( Distributions::Type a_type)

This function returns a unique integer for the Distributions::Type. For internal use when broadcasting a distribution for MPI and GPUs needs.

Parameters
a_type[in] The distribution's type.
Returns
Returns a unique integer for the distribution type.

Definition at line 238 of file MCGIDI_misc.cc.

238 {
239
240 int distributionType = 0;
241
242 switch( a_type ) {
244 distributionType = 0;
245 break;
247 distributionType = 1;
248 break;
250 distributionType = 2;
251 break;
253 distributionType = 3;
254 break;
256 distributionType = 4;
257 break;
259 distributionType = 5;
260 break;
262 distributionType = 6;
263 break;
265 distributionType = 7;
266 break;
268 distributionType = 8;
269 break;
271 distributionType = 9;
272 break;
274 distributionType = 10;
275 break;
277 distributionType = 11;
278 break;
280 distributionType = 12;
281 break;
283 distributionType = 13;
284 break;
286 distributionType = 14;
287 break;
288 }
289
290 return( distributionType );
291}

Referenced by MCGIDI::Distributions::Distribution::serialize(), and serializeDistribution().

◆ Function1dClass()

LUPI_HOST_DEVICE Function1dType MCGIDI::Function1dClass ( Functions::Function1d * funct)

Definition at line 2744 of file MCGIDI_functions.cc.

2744 {
2745
2746 if( a_function == nullptr ) return( Function1dType::none );
2747 return( a_function->type( ) );
2748}

Referenced by serializeFunction1d(), serializeFunction1d_d1(), and serializeFunction1d_d2().

◆ Function2dClass()

LUPI_HOST_DEVICE Function2dType MCGIDI::Function2dClass ( Functions::Function2d * funct)

Definition at line 3219 of file MCGIDI_functions.cc.

3219 {
3220
3221 if( a_function == nullptr ) return( Function2dType::none );
3222 return( a_function->type( ) );
3223}

Referenced by serializeFunction2d().

◆ GIDI2MCGIDI_interpolation()

LUPI_HOST_DEVICE Interpolation MCGIDI::GIDI2MCGIDI_interpolation ( ptwXY_interpolation a_interpolation)

Definition at line 2732 of file MCGIDI_functions.cc.

2732 {
2733
2734 if( a_interpolation == ptwXY_interpolationLinLin ) return( Interpolation::LINLIN );
2735 if( a_interpolation == ptwXY_interpolationLogLin ) return( Interpolation::LOGLIN );
2736 if( a_interpolation == ptwXY_interpolationLinLog ) return( Interpolation::LINLOG );
2737 if( a_interpolation == ptwXY_interpolationLogLog ) return( Interpolation::LOGLOG );
2738 if( a_interpolation == ptwXY_interpolationFlat ) return( Interpolation::FLAT );
2739 return( Interpolation::OTHER );
2740}
@ ptwXY_interpolationFlat
Definition ptwXY.h:38
@ ptwXY_interpolationLinLog
Definition ptwXY.h:37
@ ptwXY_interpolationLogLog
Definition ptwXY.h:38
@ ptwXY_interpolationLinLin
Definition ptwXY.h:37
@ ptwXY_interpolationLogLin
Definition ptwXY.h:37

Referenced by MCGIDI::Distributions::CoherentElasticTNSL::CoherentElasticTNSL(), MCGIDI::Functions::FunctionBase::FunctionBase(), MCGIDI::Functions::XYs1d::XYs1d(), and MCGIDI::Functions::XYs2d::XYs2d().

◆ GIDI_VectorDoublesToMCGIDI_VectorDoubles()

LUPI_HOST Vector< double > MCGIDI::GIDI_VectorDoublesToMCGIDI_VectorDoubles ( GIDI::Vector a_vector)
Parameters
a_vector[in] The GIDI::Vector whose contents are coped to a MCGIGI::Vector.
Returns
The MCGIGI::Vector.

Definition at line 97 of file MCGIDI_misc.cc.

97 {
98
99 Vector<double> vector( a_vector.size( ) );
100
101 for( std::size_t i1 = 0; i1 < a_vector.size( ); ++i1 ) vector[i1] = a_vector[i1];
102
103 return( vector );
104}
std::size_t size() const
Definition GIDI_data.hpp:79

Referenced by MCGIDI::MultiGroupGain::MultiGroupGain().

◆ intToDistributionType()

LUPI_HOST_DEVICE Distributions::Type MCGIDI::intToDistributionType ( int a_type)

This function returns the Distributions::Type corresponding to the integer returned by distributionTypeToInt.

Parameters
a_type[in] The value returned by distributionTypeToInt.
Returns
The Distributions::Type corresponding to a_type.

Definition at line 301 of file MCGIDI_misc.cc.

301 {
302
304
305 switch( a_type ) {
306 case 0 :
308 break;
309 case 1 :
311 break;
312 case 2 :
314 break;
315 case 3 :
317 break;
318 case 4 :
320 break;
321 case 5 :
323 break;
324 case 6 :
326 break;
327 case 7 :
329 break;
330 case 8 :
332 break;
333 case 9 :
335 break;
336 case 10 :
338 break;
339 case 11 :
341 break;
342 case 12 :
344 break;
345 case 13 :
347 break;
348 case 14 :
350 break;
351 default:
352 LUPI_THROW( "intToDistributionType: unsupported distribution type." );
353 }
354
355 return( type );
356}
#define LUPI_THROW(arg)

Referenced by MCGIDI::Distributions::Distribution::serialize(), and serializeDistribution().

◆ MCGIDI_memmove()

LUPI_HOST_DEVICE void MCGIDI::MCGIDI_memmove ( char * dest,
const char * src,
size_t n )

Definition at line 386 of file MCGIDI_string.cc.

387 {
388 // Typecast src and dest addresses to (char *)
389 char *csrc = (char *)src;
390 char *cdest = (char *)dest;
391
392 // Create a temporary array to hold data of src
393 char* temp = (char*) malloc( n);
394
395 // Copy data from csrc[] to temp[]
396 for (size_t i=0; i<n; i++)
397 temp[i] = csrc[i];
398
399 // Copy data from temp[] to cdest[]
400 for (size_t i=0; i<n; i++)
401 cdest[i] = temp[i];
402
403 free(temp);
404 }
voidp malloc(uInt size)

Referenced by MCGIDI::String::append(), MCGIDI::String::erase(), MCGIDI::String::operator+=(), MCGIDI::String::operator+=(), and MCGIDI::String::operator=().

◆ MCGIDI_popsIndex()

LUPI_HOST int MCGIDI::MCGIDI_popsIndex ( PoPI::Database const & a_pops,
std::string const & a_id )

This function returns the index in a_pops for particle a_id or -1 if a_id not in a_pops.

Parameters
a_pops[in] A PoPI::Database to retrived the particle's index from.
a_id[in] The GNDS PoPs id of the particle whose index is requested.
Returns
The index.

Definition at line 83 of file MCGIDI_misc.cc.

83 {
84
85 if( !a_pops.exists( a_id ) ) return( -1 );
86
87 return( static_cast<int>( a_pops[a_id] ) );
88}

Referenced by MCGIDI::Product::Product(), MCGIDI::Product::Product(), and MCGIDI::SetupInfo::SetupInfo().

◆ MCGIDI_popsIntid()

LUPI_HOST int MCGIDI::MCGIDI_popsIntid ( PoPI::Database const & a_pops,
std::string const & a_id )

This function returns the intid for particle a_id or -1 if a_id not in a_pops.

Parameters
a_pops[in] A PoPI::Database to retrived the particle's intid from.
a_id[in] The GNDS PoPs id of the particle whose intid is requested.
Returns
The intid.

Definition at line 51 of file MCGIDI_misc.cc.

51 {
52
55 int intid = a_pops.intid( a_id );
56
57 if( intid < 0 ) {
58 if( a_id == PoPI::IDs::neutron ) {
59 intid = PoPI::Intids::neutron; }
60 else if( a_id == PoPI::IDs::photon ) {
61 intid = PoPI::Intids::photon ; }
62 else {
63 PoPI::ParseIdInfo parseIdInfo( a_id );
64 if( parseIdInfo.isSupported( ) ) {
65 if( parseIdInfo.isNuclear( ) ) {
66 intid = 1000 * parseIdInfo.Z( ) + parseIdInfo.A( );
67 }
68 }
69 }
70 }
71 return( intid );
72}
static std::string const photon
Definition PoPI.hpp:162
static std::string const neutron
Definition PoPI.hpp:164
static std::string const FissionProductENDL99125
Definition PoPI.hpp:172
static std::string const FissionProductENDL99120
Definition PoPI.hpp:171
static int constexpr FissionProductENDL99120
Definition PoPI.hpp:180
static int constexpr photon
Definition PoPI.hpp:178
static int constexpr neutron
Definition PoPI.hpp:177
static int constexpr FissionProductENDL99125
Definition PoPI.hpp:181

Referenced by MCGIDI::OutputChannel::OutputChannel(), MCGIDI::Product::Product(), and MCGIDI::Product::Product().

◆ MCGIDI_strcmp()

LUPI_HOST_DEVICE int MCGIDI::MCGIDI_strcmp ( const char * p1,
const char * p2 )

Definition at line 357 of file MCGIDI_string.cc.

358 {
359 const unsigned char *s1 = (const unsigned char *) p1;
360 const unsigned char *s2 = (const unsigned char *) p2;
361 unsigned char c1, c2;
362
363 do
364 {
365 c1 = (unsigned char) *s1++;
366 c2 = (unsigned char) *s2++;
367 if (c1 == '\0')
368 return c1 - c2;
369 }
370 while (c1 == c2);
371
372 return c1 - c2;
373 }

Referenced by operator<(), MCGIDI::String::operator==(), and MCGIDI::String::operator==().

◆ MCGIDI_strlen()

LUPI_HOST_DEVICE size_t MCGIDI::MCGIDI_strlen ( const char * str)

Definition at line 375 of file MCGIDI_string.cc.

375 {
376 size_t len = 0;
377 while (*str != '\0') {
378 str++;
379 len++;
380 }
381 return len;
382 }

Referenced by MCGIDI::String::append(), MCGIDI::String::at(), MCGIDI::String::at(), MCGIDI::String::compare(), MCGIDI::String::operator+=(), and MCGIDI::String::operator=().

◆ MCGIDI_strncmp()

LUPI_HOST_DEVICE int MCGIDI::MCGIDI_strncmp ( const char * s1,
const char * s2,
size_t n )

Definition at line 406 of file MCGIDI_string.cc.

407 {
408 while ( n && *s1 && ( *s1 == *s2 ) )
409 {
410 ++s1;
411 ++s2;
412 --n;
413 }
414 if ( n == 0 )
415 {
416 return 0;
417 }
418 else
419 {
420 return ( *(unsigned char *)s1 - *(unsigned char *)s2 );
421 }
422 }

Referenced by MCGIDI::String::compare(), and MCGIDI::String::compare().

◆ muCOM_From_muLab()

LUPI_HOST_DEVICE int MCGIDI::muCOM_From_muLab ( double a_muLab,
double a_boostBeta,
double a_productBeta,
double & a_muPlus,
double & a_JacobianPlus,
double & a_muMinus,
double & a_JacobianMinus )

This function determines the mu value(s) in the center-of-mass frame for a specified mu value in the lab frame for a product with speed a_productBeta and a boost speed of a_productBeta. The returned value is the number of mu values in the center-of-mass frame. The return value can be 0, 1 or 2. a_muMinus and a_JacobianMinus are undefined when the returned value is less than 2. a_muPlus and a_JacobianPlus are undefined when the returned value is 0.

Parameters
a_muLab[in] The mu specified mu value in the lab frame.
a_boostBeta[in] The boost speed from the lab from to the center-of-mass frame in units of the speed-of-light.
a_productBeta[in] The speed of the product in the center-of-mass frame in units of the speed-of-light.
a_muPlus[in] The first mu value if the returned is greater than 0.
a_JacobianPlus[in] The partial derivative of mu_com with respect to mu_lab at a_muPlus.
a_muMinus[in] The second mu value if the returned value is 2.
a_JacobianMinus[in] The partial derivative of mu_com with respect to mu_lab at a_muMinus.
Returns
The number of returned center-of-mass frame mu values. Can be 0, 1 or 2.

Definition at line 185 of file MCGIDI_misc.cc.

186 {
187
188 int numberOfSolutions = 0;
189 double boostBeta2 = a_boostBeta * a_boostBeta;
190 double productBeta2 = a_productBeta * a_productBeta;
191 double muLab2 = a_muLab * a_muLab;
192 double oneMinusMuLab2 = 1.0 - muLab2;
193 double oneMinusBoostBeta2 = 1.0 - boostBeta2;
194 double oneMinusBoostBeta2MuLab2 = 1.0 - boostBeta2 * muLab2;
195
196 a_muPlus = a_muLab; // Handles case when a_productBeta is 0.0 or a_muLab is +/-1.0. Actually, when a_productBeta is 0.0 it is not defined.
197 a_muMinus = 0.0;
198
199 if( ( a_productBeta == 0.0 ) || ( a_muLab == 1.0 ) ) return( 1 );
200
201 if( a_productBeta >= a_boostBeta ) { // Intentionally treating case where a_productBeta == a_boostBeta as one solution even though is it
202 numberOfSolutions = 1; }
203 else {
204 if( a_muLab > 0.0 ) { // Only have solutions for positive mu. The next expression only test mu^2 and therefore treats negative mu like positive mu.
205 if( productBeta2 * oneMinusBoostBeta2MuLab2 > boostBeta2 * oneMinusMuLab2 ) numberOfSolutions = 2; // This ignores the case for numberOfSolutions = 1 as it probabily will never happen.
206 }
207 }
208
209 if( numberOfSolutions == 0 ) return( 0 );
210
211 double sqrt_b2minus4ac = sqrt( oneMinusBoostBeta2 * ( productBeta2 * oneMinusBoostBeta2MuLab2 - boostBeta2 * oneMinusMuLab2 ) );
212 double minusbTerm = a_boostBeta * oneMinusMuLab2;
213 double inv2a = 1.0 / ( a_productBeta * oneMinusBoostBeta2MuLab2 );
214
215 a_muPlus = ( a_muLab * sqrt_b2minus4ac - minusbTerm ) * inv2a;
216 a_muMinus = ( -a_muLab * sqrt_b2minus4ac - minusbTerm ) * inv2a; // This is meaningless when numberOfSolutions is not 1, but why add an if test.
217
218 double JacobianTerm1 = 2.0 * boostBeta2 * a_muLab / oneMinusBoostBeta2MuLab2;
219 double JacobianTerm2 = 2.0 * a_muLab * a_boostBeta / ( a_productBeta * oneMinusBoostBeta2MuLab2 );
220 double JacobianTerm3 = productBeta2 * ( 1.0 - 2.0 * boostBeta2 * muLab2 ) - boostBeta2 * ( 1.0 - 2.0 * muLab2 );
221 JacobianTerm3 *= oneMinusBoostBeta2 / ( a_productBeta * oneMinusBoostBeta2MuLab2 * sqrt_b2minus4ac );
222
223 a_JacobianPlus = fabs( a_muPlus * JacobianTerm1 + JacobianTerm2 + JacobianTerm3 );
224 a_JacobianMinus = fabs( a_muMinus * JacobianTerm1 + JacobianTerm2 - JacobianTerm3 );
225
226 return( numberOfSolutions );
227}

Referenced by MCGIDI::Distributions::AngularTwoBody::angleBiasing(), MCGIDI::Distributions::EnergyAngularMC::angleBiasing(), MCGIDI::Distributions::KalbachMann::angleBiasing(), and MCGIDI::Distributions::Uncorrelated::angleBiasing().

◆ operator+()

LUPI_HOST_DEVICE String MCGIDI::operator+ ( const String & lhs,
const String & rhs )

Definition at line 183 of file MCGIDI_string.cc.

184 {
185 return String(lhs) += rhs;
186 }

◆ operator<()

LUPI_HOST_DEVICE bool MCGIDI::operator< ( const String & s1,
const String & s2 )

Definition at line 350 of file MCGIDI_string.cc.

350 {
351 return MCGIDI_strcmp( s1.c_str(), s2.c_str() ) < 0;
352 }
LUPI_HOST_DEVICE const char * c_str() const
raw data
LUPI_HOST_DEVICE int MCGIDI_strcmp(const char *p1, const char *p2)

◆ particleKineticEnergy()

LUPI_HOST_DEVICE double MCGIDI::particleKineticEnergy ( double a_mass_unitOfEnergy,
double a_particleBeta )

This function returns a particle kinetic energy from its mass and beta (i.e., v/c) using a relativistic formula.

Parameters
a_mass_unitOfEnergy[in] The particle's mass in units of energy.
a_particleBeta[in] The particle's velocity divided by the speed of light (i.e., beta = v/c).
Returns
The relativistic kinetic energy of the particle.

Definition at line 127 of file MCGIDI_misc.cc.

127 {
128
129 if( a_particleBeta < 1e-4 ) return( 0.5 * a_mass_unitOfEnergy * a_particleBeta * a_particleBeta );
130
131 return( a_mass_unitOfEnergy * ( 1.0 / sqrt( 1.0 - a_particleBeta * a_particleBeta ) - 1.0 ) );
132}

Referenced by MCGIDI::ProtareSingle::sampleTargetBetaForUpscatterModelA().

◆ particleKineticEnergyFromBeta2()

LUPI_HOST_DEVICE double MCGIDI::particleKineticEnergyFromBeta2 ( double a_mass_unitOfEnergy,
double a_particleBeta2 )

This function is like particleKineticEnergy except that a_particleBeta2 is beta squared (i.e., (v/c)^2).

Parameters
a_mass_unitOfEnergy[in] The particle's mass in units of energy.
a_particleBeta2[in] The square of beta (i.e., beta^2 where beta = v/c).
Returns
The relativistic kinetic energy of the particle.

Definition at line 144 of file MCGIDI_misc.cc.

144 {
145
146 if( a_particleBeta2 < 1e-8 ) return( 0.5 * a_mass_unitOfEnergy * a_particleBeta2 );
147
148 return( a_mass_unitOfEnergy * ( 1.0 / sqrt( 1.0 - a_particleBeta2 ) - 1.0 ) );
149}

Referenced by MCGIDI::Distributions::AngularTwoBody::angleBiasing(), MCGIDI::Distributions::EnergyAngularMC::angleBiasing(), MCGIDI::Distributions::KalbachMann::angleBiasing(), MCGIDI::Distributions::Uncorrelated::angleBiasing(), and upScatterModelABoostParticle().

◆ ProbabilityBase1dClass()

LUPI_HOST_DEVICE ProbabilityBase1dType MCGIDI::ProbabilityBase1dClass ( Probabilities::ProbabilityBase1d * funct)

Definition at line 3288 of file MCGIDI_functions.cc.

3288 {
3289
3290 if( a_function == nullptr ) return( ProbabilityBase1dType::none );
3291 return( a_function->type( ) );
3292}

Referenced by serializeProbability1d().

◆ ProbabilityBase2dClass()

LUPI_HOST_DEVICE ProbabilityBase2dType MCGIDI::ProbabilityBase2dClass ( Probabilities::ProbabilityBase2d * funct)

Definition at line 3357 of file MCGIDI_functions.cc.

3357 {
3358
3359 if( a_function == nullptr ) return( ProbabilityBase2dType::none );
3360 return( a_function->type( ) );
3361}

Referenced by serializeProbability2d(), serializeProbability2d_d1(), and serializeProbability2d_d2().

◆ ProbabilityBase3dClass()

LUPI_HOST_DEVICE ProbabilityBase3dType MCGIDI::ProbabilityBase3dClass ( Probabilities::ProbabilityBase3d * funct)

Definition at line 4063 of file MCGIDI_functions.cc.

4063 {
4064
4065 if( a_function == nullptr ) return( ProbabilityBase3dType::none );
4066 return( a_function->type( ) );
4067}

Referenced by serializeProbability3d().

◆ protareFromGIDIProtare()

LUPI_HOST Protare * MCGIDI::protareFromGIDIProtare ( LUPI::StatusMessageReporting & a_smr,
GIDI::Protare 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,
bool a_allowFixedGrid )

Returns the proper MCGIDI protare base on the type of GIDI protare.

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

33 {
34
35 Protare *protare( nullptr );
36
37 if( a_protare.protareType( ) == GIDI::ProtareType::single ) {
38 protare = new ProtareSingle( a_smr, static_cast<GIDI::ProtareSingle const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
39 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, a_allowFixedGrid ); }
40 else if( a_protare.protareType( ) == GIDI::ProtareType::composite ) {
41 protare = new ProtareComposite( a_smr, static_cast<GIDI::ProtareComposite const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
42 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, false ); }
43 else if( a_protare.protareType( ) == GIDI::ProtareType::TNSL ) {
44 protare = new ProtareTNSL( a_smr, static_cast<GIDI::ProtareTNSL const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
45 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, false );
46 }
47
48 return( protare );
49}

Referenced by G4GIDI::readTarget().

◆ serializeACE_URR_probabilityTables()

LUPI_HOST_DEVICE ACE_URR_probabilityTables * MCGIDI::serializeACE_URR_probabilityTables ( ACE_URR_probabilityTables * a_ACE_URR_probabilityTables,
LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )

This method serializes a ACE_URR_probabilityTables instance pointed to by a_ACE_URR_probabilityTables.

Parameters
a_ACE_URR_probabilityTables[in] Specifies the action of this method.
a_buffer[in] The buffer to read or write data to depending on a_mode.
a_mode[in] Specifies the action of this method.
Returns
A pointer to the serialized ACE_URR_probabilityTables instance.

Definition at line 379 of file MCGIDI_URR.cc.

380 {
381
382 int type = 0;
383 if( a_ACE_URR_probabilityTables != nullptr ) type = 1;
384 DATA_MEMBER_INT( type, a_buffer, a_mode );
385 if( type == 0 ) return( nullptr );
386
387 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
388 if( a_buffer.m_placement != nullptr ) {
389 a_ACE_URR_probabilityTables = new(a_buffer.m_placement) ACE_URR_probabilityTables;
390 a_buffer.incrementPlacement( sizeof( ACE_URR_probabilityTables ) ); }
391 else {
392 a_ACE_URR_probabilityTables = new ACE_URR_probabilityTables;
393 }
394 }
395
397
398 a_ACE_URR_probabilityTables->serialize( a_buffer, a_mode );
399
400 return( a_ACE_URR_probabilityTables );
401}
#define DATA_MEMBER_INT( member, buf, mode)
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

Referenced by MCGIDI::HeatedCrossSectionContinuousEnergy::serialize(), and MCGIDI::HeatedReactionCrossSectionContinuousEnergy::serialize().

◆ serializeDelayedNeutrons()

LUPI_HOST_DEVICE void MCGIDI::serializeDelayedNeutrons ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Vector< DelayedNeutron * > & a_delayedNeutrons )

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

Parameters
a_delayedNeutrons[in] The delayed neutrons to serialize.
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 404 of file MCGIDI_misc.cc.

404 {
405
406 std::size_t vectorSize = a_delayedNeutrons.size( );
407 int vectorSizeInt = (int) vectorSize;
408 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
409 vectorSize = (std::size_t) vectorSizeInt;
410
411 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
412 a_delayedNeutrons.resize( vectorSize, &a_buffer.m_placement );
413 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
414 if( a_buffer.m_placement != nullptr ) {
415 a_delayedNeutrons[vectorIndex] = new(a_buffer.m_placement) DelayedNeutron;
416 a_buffer.incrementPlacement( sizeof( DelayedNeutron ) );
417 }
418 else {
419 a_delayedNeutrons[vectorIndex] = new DelayedNeutron;
420 }
421 } }
422 else if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
423 a_buffer.m_placement += a_delayedNeutrons.internalSize( );
424 a_buffer.incrementPlacement( sizeof( DelayedNeutron ) * vectorSize );
425 }
426
427 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
428 a_delayedNeutrons[vectorIndex]->serialize( a_buffer, a_mode );
429 }
430}
LUPI_HOST_DEVICE std::size_t internalSize() const
LUPI_HOST_DEVICE void resize(std::size_t s, char **address=nullptr, bool mem_flag=CPU_MEM)

Referenced by MCGIDI::OutputChannel::serialize(), and MCGIDI::Reaction::serialize().

◆ serializeDistribution()

LUPI_HOST_DEVICE Distributions::Distribution * MCGIDI::serializeDistribution ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Distributions::Distribution * a_distribution )

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
a_mode[in] Specifies the action of this method.

Definition at line 1535 of file MCGIDI_distributions.cc.

1536 {
1537
1539 if( a_distribution != nullptr ) type = a_distribution->type( );
1540 int distributionType = distributionTypeToInt( type );
1541 DATA_MEMBER_INT( distributionType, a_buffer, a_mode );
1542 type = intToDistributionType( distributionType );
1543
1544 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
1545 switch( type ) {
1547 a_distribution = nullptr;
1548 break;
1550 if (a_buffer.m_placement != nullptr) {
1551 a_distribution = new(a_buffer.m_placement) Distributions::Unspecified;
1552 a_buffer.incrementPlacement( sizeof( Distributions::Unspecified ) ); }
1553 else {
1554 a_distribution = new Distributions::Unspecified;
1555 }
1556 break;
1558 if (a_buffer.m_placement != nullptr) {
1559 a_distribution = new(a_buffer.m_placement) Distributions::AngularTwoBody;
1560 a_buffer.incrementPlacement( sizeof( Distributions::AngularTwoBody ) ); }
1561 else {
1562 a_distribution = new Distributions::AngularTwoBody;
1563 }
1564 break;
1566 if (a_buffer.m_placement != nullptr) {
1567 a_distribution = new(a_buffer.m_placement) Distributions::KalbachMann;
1568 a_buffer.incrementPlacement( sizeof( Distributions::KalbachMann ) ); }
1569 else {
1570 a_distribution = new Distributions::KalbachMann;
1571 }
1572 break;
1574 if (a_buffer.m_placement != nullptr) {
1575 a_distribution = new(a_buffer.m_placement) Distributions::Uncorrelated;
1576 a_buffer.incrementPlacement( sizeof( Distributions::Uncorrelated ) ); }
1577 else {
1578 a_distribution = new Distributions::Uncorrelated;
1579 }
1580 break;
1582 if (a_buffer.m_placement != nullptr) {
1583 a_distribution = new(a_buffer.m_placement) Distributions::Branching3d;
1584 a_buffer.incrementPlacement( sizeof( Distributions::Branching3d ) ); }
1585 else {
1586 a_distribution = new Distributions::Branching3d;
1587 }
1588 break;
1590 if (a_buffer.m_placement != nullptr) {
1591 a_distribution = new(a_buffer.m_placement) Distributions::EnergyAngularMC;
1592 a_buffer.incrementPlacement( sizeof( Distributions::EnergyAngularMC ) ); }
1593 else {
1594 a_distribution = new Distributions::EnergyAngularMC;
1595 }
1596 break;
1598 if (a_buffer.m_placement != nullptr) {
1599 a_distribution = new(a_buffer.m_placement) Distributions::AngularEnergyMC;
1600 a_buffer.incrementPlacement( sizeof( Distributions::AngularEnergyMC ) ); }
1601 else {
1602 a_distribution = new Distributions::AngularEnergyMC;
1603 }
1604 break;
1606 if( a_buffer.m_placement != nullptr ) {
1607 a_distribution = new(a_buffer.m_placement) Distributions::CoherentPhotoAtomicScattering;
1609 else {
1611 }
1612 break;
1614 if( a_buffer.m_placement != nullptr ) {
1615 a_distribution = new(a_buffer.m_placement) Distributions::IncoherentPhotoAtomicScattering;
1617 else {
1619 }
1620 break;
1622 if( a_buffer.m_placement != nullptr ) {
1625 else {
1627 }
1628 break;
1630 if( a_buffer.m_placement != nullptr ) {
1631 a_distribution = new(a_buffer.m_placement) Distributions::PairProductionGamma;
1633 else {
1634 a_distribution = new Distributions::PairProductionGamma;
1635 }
1636 break;
1638 if( a_buffer.m_placement != nullptr ) {
1639 a_distribution = new(a_buffer.m_placement) Distributions::CoherentElasticTNSL;
1641 else {
1642 a_distribution = new Distributions::CoherentElasticTNSL;
1643 }
1644 break;
1646 if( a_buffer.m_placement != nullptr ) {
1647 a_distribution = new(a_buffer.m_placement) Distributions::IncoherentElasticTNSL;
1649 else {
1650 a_distribution = new Distributions::IncoherentElasticTNSL;
1651 }
1652 break;
1654 if( a_buffer.m_placement != nullptr ) {
1657 else {
1659 }
1660 break;
1661 }
1662 }
1663
1664 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
1665 switch( type ) {
1667 break;
1669 a_buffer.incrementPlacement( sizeof( Distributions::Unspecified ) );
1670 break;
1673 break;
1675 a_buffer.incrementPlacement( sizeof( Distributions::KalbachMann ) );
1676 break;
1678 a_buffer.incrementPlacement( sizeof( Distributions::Uncorrelated ) );
1679 break;
1681 a_buffer.incrementPlacement( sizeof( Distributions::Branching3d ) );
1682 break;
1685 break;
1688 break;
1691 break;
1694 break;
1697 break;
1700 break;
1703 break;
1706 break;
1709 break;
1710 }
1711 }
1712
1713 switch( type ) {
1715 break;
1717 static_cast<Distributions::Unspecified *>( a_distribution )->serialize( a_buffer, a_mode );
1718 break;
1720 static_cast<Distributions::AngularTwoBody *>( a_distribution )->serialize( a_buffer, a_mode );
1721 break;
1723 static_cast<Distributions::KalbachMann *>( a_distribution )->serialize( a_buffer, a_mode );
1724 break;
1726 static_cast<Distributions::Uncorrelated *>( a_distribution )->serialize( a_buffer, a_mode );
1727 break;
1729 static_cast<Distributions::Branching3d *>( a_distribution )->serialize( a_buffer, a_mode );
1730 break;
1732 static_cast<Distributions::EnergyAngularMC *>( a_distribution )->serialize( a_buffer, a_mode );
1733 break;
1735 static_cast<Distributions::AngularEnergyMC *>( a_distribution )->serialize( a_buffer, a_mode );
1736 break;
1738 static_cast<Distributions::CoherentPhotoAtomicScattering *>( a_distribution )->serialize( a_buffer, a_mode );
1739 break;
1741 static_cast<Distributions::IncoherentPhotoAtomicScattering *>( a_distribution )->serialize( a_buffer, a_mode );
1742 break;
1744 static_cast<Distributions::IncoherentBoundToFreePhotoAtomicScattering *>( a_distribution )->serialize( a_buffer, a_mode );
1745 break;
1747 static_cast<Distributions::PairProductionGamma *>( a_distribution )->serialize( a_buffer, a_mode );
1748 break;
1750 static_cast<Distributions::CoherentElasticTNSL *>( a_distribution )->serialize( a_buffer, a_mode );
1751 break;
1753 static_cast<Distributions::IncoherentElasticTNSL *>( a_distribution )->serialize( a_buffer, a_mode );
1754 break;
1756 static_cast<Distributions::IncoherentPhotoAtomicScatteringElectron *>( a_distribution )->serialize( a_buffer, a_mode );
1757 break;
1758 }
1759
1760 return( a_distribution );
1761}
LUPI_HOST_DEVICE int distributionTypeToInt(Distributions::Type a_type)
LUPI_HOST_DEVICE Distributions::Type intToDistributionType(int a_type)

Referenced by MCGIDI::Product::serialize().

◆ serializeFissionResiduals()

LUPI_HOST_DEVICE void MCGIDI::serializeFissionResiduals ( GIDI::Construction::FissionResiduals & a_fissionResiduals,
LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )
Parameters
a_fissionResiduals[in] A reference to the GIDI::Construction::FissionResiduals reference serialize.
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 467 of file MCGIDI_misc.cc.

468 {
469
470 int fissionResidualsInt = 0;
471
472 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
473 switch( a_fissionResiduals ) {
475 break;
477 fissionResidualsInt = 1;
478 break;
480 fissionResidualsInt = 2;
481 break;
482 }
483 }
484
485 DATA_MEMBER_INT( fissionResidualsInt, a_buffer, a_mode );
486
487 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
488 switch( fissionResidualsInt ) {
489 case 0 :
491 break;
492 case 1 :
494 break;
495 case 2 :
497 break;
498 }
499 }
500}

Referenced by MCGIDI::Reaction::serialize().

◆ serializeFunction1d()

LUPI_HOST_DEVICE Functions::Function1d * MCGIDI::serializeFunction1d ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Functions::Function1d * a_function1d )

Definition at line 2753 of file MCGIDI_functions.cc.

2754 {
2755
2756 int type = 0;
2757
2758 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
2759 Function1dType fType = Function1dClass( a_function1d );
2760
2761 switch( fType ) {
2763 break;
2765 type = 1;
2766 break;
2767 case Function1dType::XYs :
2768 type = 2;
2769 break;
2771 type = 3;
2772 break;
2774 type = 4;
2775 break;
2777 type = 5;
2778 break;
2780 type = 6;
2781 break;
2783 type = 7;
2784 break;
2785 default:
2786 String message( "serializeFunction1d: Unsupported Function1d: " + a_function1d->typeString( ) );
2787 LUPI_THROW( message.c_str( ) );
2788 }
2789 }
2790
2791 DATA_MEMBER_INT( type, a_buffer, a_mode );
2792
2793 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2794 a_function1d = nullptr;
2795 switch( type ) {
2796 case 0 :
2797 break;
2798 case 1 :
2799 if( a_buffer.m_placement != nullptr ) {
2800 a_function1d = new(a_buffer.m_placement) Functions::Constant1d;
2801 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) ); }
2802 else {
2803 a_function1d = new Functions::Constant1d;
2804 }
2805 break;
2806 case 2 :
2807 if( a_buffer.m_placement != nullptr ) {
2808 a_function1d = new(a_buffer.m_placement) Functions::XYs1d;
2809 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) ); }
2810 else {
2811 a_function1d = new Functions::XYs1d;
2812 }
2813 break;
2814 case 3 :
2815 if( a_buffer.m_placement != nullptr ) {
2816 a_function1d = new(a_buffer.m_placement) Functions::Polynomial1d;
2817 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) ); }
2818 else {
2819 a_function1d = new Functions::Polynomial1d;
2820 }
2821 break;
2822 case 4 :
2823 if( a_buffer.m_placement != nullptr ) {
2824 a_function1d = new(a_buffer.m_placement) Functions::Gridded1d;
2825 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) ); }
2826 else {
2827 a_function1d = new Functions::Gridded1d;
2828 }
2829 break;
2830 case 5 :
2831 if( a_buffer.m_placement != nullptr ) {
2832 a_function1d = new(a_buffer.m_placement) Functions::Regions1d;
2833 a_buffer.incrementPlacement( sizeof( Functions::Regions1d ) ); }
2834 else {
2835 a_function1d = new Functions::Regions1d;
2836 }
2837 break;
2838 case 6 :
2839 if( a_buffer.m_placement != nullptr ) {
2840 a_function1d = new(a_buffer.m_placement) Functions::Branching1d;
2841 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) ); }
2842 else {
2843 a_function1d = new Functions::Branching1d;
2844 }
2845 break;
2846 case 7 :
2847 if( a_buffer.m_placement != nullptr ) {
2850 else {
2852 }
2853 break;
2854 default: // This should never happen as Unpack should be called after Pack which checks type.
2855 LUPI_THROW( "serializeFunction1d: Unsupported Function1d:" );
2856 }
2857 }
2858
2859 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2860 switch( type ) {
2861 case 0 :
2862 break;
2863 case 1 :
2864 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) );
2865 break;
2866 case 2 :
2867 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) );
2868 break;
2869 case 3 :
2870 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) );
2871 break;
2872 case 4 :
2873 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) );
2874 break;
2875 case 5 :
2876 a_buffer.incrementPlacement( sizeof( Functions::Regions1d ) );
2877 break;
2878 case 6 :
2879 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) );
2880 break;
2881 case 7 :
2883 break;
2884 default:
2885 LUPI_THROW( "serializeFunction1d: Unsupported Function1d:" );
2886 }
2887 }
2888
2889 if( a_function1d != nullptr ) {
2890 switch( type ) {
2891 case 0 :
2892 break;
2893 case 1 :
2894 static_cast<Functions::Constant1d *>( a_function1d )->serialize( a_buffer, a_mode );
2895 break;
2896 case 2 :
2897 static_cast<Functions::XYs1d *>( a_function1d )->serialize( a_buffer, a_mode );
2898 break;
2899 case 3 :
2900 static_cast<Functions::Polynomial1d *>( a_function1d )->serialize( a_buffer, a_mode );
2901 break;
2902 case 4 :
2903 static_cast<Functions::Gridded1d *>( a_function1d )->serialize( a_buffer, a_mode );
2904 break;
2905 case 5 :
2906 static_cast<Functions::Regions1d *>( a_function1d )->serialize( a_buffer, a_mode );
2907 break;
2908 case 6 :
2909 static_cast<Functions::Branching1d *>( a_function1d )->serialize( a_buffer, a_mode );
2910 break;
2911 case 7 :
2912 static_cast<Functions::TerrellFissionNeutronMultiplicityModel *>( a_function1d )->serialize( a_buffer, a_mode );
2913 break;
2914 default:
2915 LUPI_THROW( "serializeFunction1d: Unsupported Function1d:" );
2916 }
2917 }
2918
2919 return( a_function1d );
2920}
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE Function1dType Function1dClass(Functions::Function1d *funct)

Referenced by MCGIDI::OutputChannel::serialize(), MCGIDI::Product::serialize(), and MCGIDI::Reaction::serialize().

◆ serializeFunction1d_d1()

LUPI_HOST_DEVICE Functions::Function1d_d1 * MCGIDI::serializeFunction1d_d1 ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Functions::Function1d_d1 * a_function1d )

Definition at line 2925 of file MCGIDI_functions.cc.

2926 {
2927
2928 int type = 0;
2929
2930 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
2931 Function1dType fType = Function1dClass( a_function1d );
2932
2933 switch( fType ) {
2935 break;
2937 type = 1;
2938 break;
2939 case Function1dType::XYs :
2940 type = 2;
2941 break;
2943 type = 3;
2944 break;
2946 type = 4;
2947 break;
2949 type = 5;
2950 break;
2952 type = 6;
2953 break;
2954 default:
2955 String message( "serializeFunction1d_d1: Unsupported Function1d: " + a_function1d->typeString( ) );
2956 LUPI_THROW( message.c_str( ) );
2957 }
2958 }
2959
2960 DATA_MEMBER_INT( type, a_buffer, a_mode );
2961
2962 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2963 a_function1d = nullptr;
2964 switch( type ) {
2965 case 0 :
2966 break;
2967 case 1 :
2968 if( a_buffer.m_placement != nullptr ) {
2969 a_function1d = new(a_buffer.m_placement) Functions::Constant1d;
2970 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) ); }
2971 else {
2972 a_function1d = new Functions::Constant1d;
2973 }
2974 break;
2975 case 2 :
2976 if( a_buffer.m_placement != nullptr ) {
2977 a_function1d = new(a_buffer.m_placement) Functions::XYs1d;
2978 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) ); }
2979 else {
2980 a_function1d = new Functions::XYs1d;
2981 }
2982 break;
2983 case 3 :
2984 if( a_buffer.m_placement != nullptr ) {
2985 a_function1d = new(a_buffer.m_placement) Functions::Polynomial1d;
2986 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) ); }
2987 else {
2988 a_function1d = new Functions::Polynomial1d;
2989 }
2990 break;
2991 case 4 :
2992 if( a_buffer.m_placement != nullptr ) {
2993 a_function1d = new(a_buffer.m_placement) Functions::Gridded1d;
2994 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) ); }
2995 else {
2996 a_function1d = new Functions::Gridded1d;
2997 }
2998 break;
2999 case 5 :
3000 if( a_buffer.m_placement != nullptr ) {
3001 a_function1d = new(a_buffer.m_placement) Functions::Regions1d;
3002 a_buffer.incrementPlacement( sizeof( Functions::Regions1d ) ); }
3003 else {
3004 a_function1d = new Functions::Regions1d;
3005 }
3006 break;
3007 case 6 :
3008 if( a_buffer.m_placement != nullptr ) {
3009 a_function1d = new(a_buffer.m_placement) Functions::Branching1d;
3010 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) ); }
3011 else {
3012 a_function1d = new Functions::Branching1d;
3013 }
3014 break;
3015 default: // This should never happen as Unpack should be called after Pack which checks type.
3016 LUPI_THROW( "serializeFunction1d_d1: Unsupported Function1d:" );
3017 }
3018 }
3019
3020 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3021 switch( type ) {
3022 case 0 :
3023 break;
3024 case 1 :
3025 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) );
3026 break;
3027 case 2 :
3028 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) );
3029 break;
3030 case 3 :
3031 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) );
3032 break;
3033 case 4 :
3034 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) );
3035 break;
3036 case 5 :
3037 a_buffer.incrementPlacement( sizeof( Functions::Regions1d ) );
3038 break;
3039 case 6 :
3040 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) );
3041 break;
3042 default:
3043 LUPI_THROW( "serializeFunction1d_d1: Unsupported Function1d:" );
3044 }
3045 }
3046
3047 if( a_function1d != nullptr ) {
3048 switch( type ) {
3049 case 0 :
3050 break;
3051 case 1 :
3052 static_cast<Functions::Constant1d *>( a_function1d )->serialize( a_buffer, a_mode );
3053 break;
3054 case 2 :
3055 static_cast<Functions::XYs1d *>( a_function1d )->serialize( a_buffer, a_mode );
3056 break;
3057 case 3 :
3058 static_cast<Functions::Polynomial1d *>( a_function1d )->serialize( a_buffer, a_mode );
3059 break;
3060 case 4 :
3061 static_cast<Functions::Gridded1d *>( a_function1d )->serialize( a_buffer, a_mode );
3062 break;
3063 case 5 :
3064 static_cast<Functions::Regions1d *>( a_function1d )->serialize( a_buffer, a_mode );
3065 break;
3066 case 6 :
3067 static_cast<Functions::Branching1d *>( a_function1d )->serialize( a_buffer, a_mode );
3068 break;
3069 default:
3070 LUPI_THROW( "serializeFunction1d_d1: Unsupported Function1d:" );
3071 }
3072 }
3073
3074 return( a_function1d );
3075}

Referenced by MCGIDI::Distributions::CoherentPhotoAtomicScattering::serialize(), MCGIDI::Distributions::IncoherentElasticTNSL::serialize(), MCGIDI::Functions::TerrellFissionNeutronMultiplicityModel::serialize(), MCGIDI::Functions::XYs2d::serialize(), MCGIDI::OutputChannel::serialize(), MCGIDI::Probabilities::Evaporation2d::serialize(), MCGIDI::Probabilities::GeneralEvaporation2d::serialize(), MCGIDI::Probabilities::SimpleMaxwellianFission2d::serialize(), MCGIDI::Probabilities::Watt2d::serialize(), MCGIDI::Probabilities::WeightedFunctionals2d::serialize(), and serializeQs().

◆ serializeFunction1d_d2()

LUPI_HOST_DEVICE Functions::Function1d_d2 * MCGIDI::serializeFunction1d_d2 ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Functions::Function1d_d2 * a_function1d )

Definition at line 3080 of file MCGIDI_functions.cc.

3081 {
3082
3083 int type = 0;
3084
3085 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3086 Function1dType fType = Function1dClass( a_function1d );
3087
3088 switch( fType ) {
3090 break;
3092 type = 1;
3093 break;
3094 case Function1dType::XYs :
3095 type = 2;
3096 break;
3098 type = 3;
3099 break;
3101 type = 4;
3102 break;
3103 break;
3105 type = 6;
3106 break;
3107 default:
3108 String message( "serializeFunction1d_d2: Unsupported Function1d: " + a_function1d->typeString( ) );
3109 LUPI_THROW( message.c_str( ) );
3110 }
3111 }
3112
3113 DATA_MEMBER_INT( type, a_buffer, a_mode );
3114
3115 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3116 a_function1d = nullptr;
3117 switch( type ) {
3118 case 0 :
3119 break;
3120 case 1 :
3121 if( a_buffer.m_placement != nullptr ) {
3122 a_function1d = new(a_buffer.m_placement) Functions::Constant1d;
3123 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) ); }
3124 else {
3125 a_function1d = new Functions::Constant1d;
3126 }
3127 break;
3128 case 2 :
3129 if( a_buffer.m_placement != nullptr ) {
3130 a_function1d = new(a_buffer.m_placement) Functions::XYs1d;
3131 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) ); }
3132 else {
3133 a_function1d = new Functions::XYs1d;
3134 }
3135 break;
3136 case 3 :
3137 if( a_buffer.m_placement != nullptr ) {
3138 a_function1d = new(a_buffer.m_placement) Functions::Polynomial1d;
3139 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) ); }
3140 else {
3141 a_function1d = new Functions::Polynomial1d;
3142 }
3143 break;
3144 case 4 :
3145 if( a_buffer.m_placement != nullptr ) {
3146 a_function1d = new(a_buffer.m_placement) Functions::Gridded1d;
3147 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) ); }
3148 else {
3149 a_function1d = new Functions::Gridded1d;
3150 }
3151 break;
3152 case 6 :
3153 if( a_buffer.m_placement != nullptr ) {
3154 a_function1d = new(a_buffer.m_placement) Functions::Branching1d;
3155 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) ); }
3156 else {
3157 a_function1d = new Functions::Branching1d;
3158 }
3159 break;
3160 default: // This should never happen as Unpack should be called after Pack which checks type.
3161 LUPI_THROW( "serializeFunction1d_d2: Unsupported Function1d:" );
3162 }
3163 }
3164
3165 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3166 switch( type ) {
3167 case 0 :
3168 break;
3169 case 1 :
3170 a_buffer.incrementPlacement( sizeof( Functions::Constant1d ) );
3171 break;
3172 case 2 :
3173 a_buffer.incrementPlacement( sizeof( Functions::XYs1d ) );
3174 break;
3175 case 3 :
3176 a_buffer.incrementPlacement( sizeof( Functions::Polynomial1d ) );
3177 break;
3178 case 4 :
3179 a_buffer.incrementPlacement( sizeof( Functions::Gridded1d ) );
3180 break;
3181 case 6 :
3182 a_buffer.incrementPlacement( sizeof( Functions::Branching1d ) );
3183 break;
3184 default:
3185 LUPI_THROW( "serializeFunction1d_d2: Unsupported Function1d:" );
3186 }
3187 }
3188
3189 if( a_function1d != nullptr ) {
3190 switch( type ) {
3191 case 0 :
3192 break;
3193 case 1 :
3194 static_cast<Functions::Constant1d *>( a_function1d )->serialize( a_buffer, a_mode );
3195 break;
3196 case 2 :
3197 static_cast<Functions::XYs1d *>( a_function1d )->serialize( a_buffer, a_mode );
3198 break;
3199 case 3 :
3200 static_cast<Functions::Polynomial1d *>( a_function1d )->serialize( a_buffer, a_mode );
3201 break;
3202 case 4 :
3203 static_cast<Functions::Gridded1d *>( a_function1d )->serialize( a_buffer, a_mode );
3204 break;
3205 case 6 :
3206 static_cast<Functions::Branching1d *>( a_function1d )->serialize( a_buffer, a_mode );
3207 break;
3208 default:
3209 LUPI_THROW( "serializeFunction1d_d2: Unsupported Function1d:" );
3210 }
3211 }
3212
3213 return( a_function1d );
3214}

Referenced by MCGIDI::Functions::Regions1d::serialize().

◆ serializeFunction2d()

LUPI_HOST_DEVICE Functions::Function2d * MCGIDI::serializeFunction2d ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Functions::Function2d * a_function2d )

Definition at line 3228 of file MCGIDI_functions.cc.

3228 {
3229
3230 int type = 0;
3231
3232 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3233 Function2dType fType = Function2dClass( a_function2d );
3234
3235 switch( fType ) {
3237 break;
3238 case Function2dType::XYs :
3239 type = 1;
3240 break;
3241 }
3242 }
3243
3244 DATA_MEMBER_INT( type, a_buffer, a_mode );
3245
3246 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3247 switch( type ) {
3248 case 0 :
3249 a_function2d = nullptr;
3250 break;
3251 case 1 :
3252 if( a_buffer.m_placement != nullptr ) {
3253 a_function2d = new(a_buffer.m_placement) Functions::XYs2d;
3254 a_buffer.incrementPlacement( sizeof( Functions::XYs2d ) ); }
3255 else {
3256 a_function2d = new Functions::XYs2d;
3257 }
3258 break;
3259 }
3260 }
3261
3262 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3263 switch( type ) {
3264 case 0 :
3265 break;
3266 case 1 :
3267 a_buffer.incrementPlacement( sizeof( Functions::XYs2d ) );
3268 break;
3269 }
3270 }
3271
3272 if( a_function2d != nullptr ) {
3273 switch( type ) {
3274 case 0 :
3275 break;
3276 case 1 :
3277 static_cast<Functions::XYs2d *>( a_function2d )->serialize( a_buffer, a_mode );
3278 break;
3279 }
3280 }
3281
3282 return( a_function2d );
3283}
LUPI_HOST_DEVICE Function2dType Function2dClass(Functions::Function2d *funct)

Referenced by MCGIDI::Distributions::KalbachMann::serialize().

◆ serializeProbability1d()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase1d * MCGIDI::serializeProbability1d ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Probabilities::ProbabilityBase1d * a_probability1d )

Definition at line 3297 of file MCGIDI_functions.cc.

3298 {
3299
3300 int type = 0;
3301
3302 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3303 ProbabilityBase1dType pType = ProbabilityBase1dClass( a_probability1d );
3304
3305 switch( pType ) {
3307 break;
3309 type = 1;
3310 break;
3311 }
3312 }
3313
3314 DATA_MEMBER_INT( type, a_buffer, a_mode );
3315
3316 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3317 switch( type ) {
3318 case 0 :
3319 a_probability1d = nullptr;
3320 break;
3321 case 1 :
3322 if( a_buffer.m_placement != nullptr ) {
3323 a_probability1d = new(a_buffer.m_placement) Probabilities::Xs_pdf_cdf1d;
3324 a_buffer.incrementPlacement( sizeof( Probabilities::Xs_pdf_cdf1d ) ); }
3325 else {
3326 a_probability1d = new Probabilities::Xs_pdf_cdf1d;
3327 }
3328 break;
3329 }
3330 }
3331 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3332 switch( type ) {
3333 case 0 :
3334 break;
3335 case 1 :
3336 a_buffer.incrementPlacement( sizeof( Probabilities::Xs_pdf_cdf1d ) );
3337 break;
3338 }
3339 }
3340
3341 if( a_probability1d != nullptr ) {
3342 switch( type ) {
3343 case 0 :
3344 break;
3345 case 1 :
3346 static_cast<Probabilities::Xs_pdf_cdf1d *>( a_probability1d )->serialize( a_buffer, a_mode );
3347 break;
3348 }
3349 }
3350
3351 return( a_probability1d );
3352}
LUPI_HOST_DEVICE ProbabilityBase1dType ProbabilityBase1dClass(Probabilities::ProbabilityBase1d *funct)

Referenced by MCGIDI::Probabilities::GeneralEvaporation2d::serialize(), MCGIDI::Probabilities::NBodyPhaseSpace2d::serialize(), and MCGIDI::Probabilities::XYs2d::serialize().

◆ serializeProbability2d()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d * MCGIDI::serializeProbability2d ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Probabilities::ProbabilityBase2d * a_probability2d )

Definition at line 3366 of file MCGIDI_functions.cc.

3367 {
3368
3369 int type = 0;
3370
3371 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3372 ProbabilityBase2dType pType = ProbabilityBase2dClass( a_probability2d );
3373
3374 switch( pType ) {
3376 break;
3378 type = 1;
3379 break;
3381 type = 2;
3382 break;
3384 type = 3;
3385 break;
3387 type = 4;
3388 break;
3390 type = 5;
3391 break;
3393 type = 6;
3394 break;
3396 type = 7;
3397 break;
3399 type = 8;
3400 break;
3402 type = 9;
3403 break;
3405 type = 10;
3406 break;
3408 type = 11;
3409 break;
3411 type = 12;
3412 break;
3413 }
3414 }
3415
3416 DATA_MEMBER_INT( type, a_buffer, a_mode );
3417
3418 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3419 switch( type ) {
3420 case 0 :
3421 a_probability2d = nullptr;
3422 break;
3423 case 1 :
3424 if( a_buffer.m_placement != nullptr ) {
3425 a_probability2d = new(a_buffer.m_placement) Probabilities::XYs2d;
3426 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) ); }
3427 else {
3428 a_probability2d = new Probabilities::XYs2d;
3429 }
3430 break;
3431 case 2 :
3432 if( a_buffer.m_placement != nullptr ) {
3433 a_probability2d = new(a_buffer.m_placement) Probabilities::Regions2d;
3434 a_buffer.incrementPlacement( sizeof( Probabilities::Regions2d ) ); }
3435 else {
3436 a_probability2d = new Probabilities::Regions2d;
3437 }
3438 break;
3439 case 3 :
3440 if( a_buffer.m_placement != nullptr ) {
3441 a_probability2d = new(a_buffer.m_placement) Probabilities::Isotropic2d;
3442 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) ); }
3443 else {
3444 a_probability2d = new Probabilities::Isotropic2d;
3445 }
3446 break;
3447 case 4 :
3448 if( a_buffer.m_placement != nullptr ) {
3449 a_probability2d = new(a_buffer.m_placement) Probabilities::DiscreteGamma2d;
3450 a_buffer.incrementPlacement( sizeof( Probabilities::DiscreteGamma2d ) ); }
3451 else {
3452 a_probability2d = new Probabilities::DiscreteGamma2d;
3453 }
3454 break;
3455 case 5 :
3456 if( a_buffer.m_placement != nullptr ) {
3457 a_probability2d = new(a_buffer.m_placement) Probabilities::PrimaryGamma2d;
3458 a_buffer.incrementPlacement( sizeof( Probabilities::PrimaryGamma2d ) ); }
3459 else {
3460 a_probability2d = new Probabilities::PrimaryGamma2d;
3461 }
3462 break;
3463 case 6 :
3464 if( a_buffer.m_placement != nullptr ) {
3465 a_probability2d = new(a_buffer.m_placement) Probabilities::Recoil2d;
3466 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) ); }
3467 else {
3468 a_probability2d = new Probabilities::Recoil2d;
3469 }
3470 break;
3471 case 7 :
3472 if( a_buffer.m_placement != nullptr ) {
3473 a_probability2d = new(a_buffer.m_placement) Probabilities::NBodyPhaseSpace2d;
3475 else {
3476 a_probability2d = new Probabilities::NBodyPhaseSpace2d;
3477 }
3478 break;
3479 case 8 :
3480 if( a_buffer.m_placement != nullptr ) {
3481 a_probability2d = new(a_buffer.m_placement) Probabilities::Evaporation2d;
3482 a_buffer.incrementPlacement( sizeof( Probabilities::Evaporation2d ) ); }
3483 else {
3484 a_probability2d = new Probabilities::Evaporation2d;
3485 }
3486 break;
3487 case 9 :
3488 if( a_buffer.m_placement != nullptr ) {
3489 a_probability2d = new(a_buffer.m_placement) Probabilities::GeneralEvaporation2d;
3491 else {
3492 a_probability2d = new Probabilities::GeneralEvaporation2d;
3493 }
3494 break;
3495 case 10 :
3496 if( a_buffer.m_placement != nullptr ) {
3497 a_probability2d = new(a_buffer.m_placement) Probabilities::SimpleMaxwellianFission2d;
3499 else {
3500 a_probability2d = new Probabilities::SimpleMaxwellianFission2d;
3501 }
3502 break;
3503 case 11 :
3504 if( a_buffer.m_placement != nullptr ) {
3505 a_probability2d = new(a_buffer.m_placement) Probabilities::Watt2d;
3506 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) ); }
3507 else {
3508 a_probability2d = new Probabilities::Watt2d;
3509 }
3510 break;
3511 case 12 :
3512 if( a_buffer.m_placement != nullptr ) {
3513 a_probability2d = new(a_buffer.m_placement) Probabilities::WeightedFunctionals2d;
3515 else {
3516 a_probability2d = new Probabilities::WeightedFunctionals2d;
3517 }
3518 break;
3519 }
3520 }
3521
3522 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3523 switch( type ) {
3524 case 0 :
3525 break;
3526 case 1 :
3527 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) );
3528 break;
3529 case 2 :
3530 a_buffer.incrementPlacement( sizeof( Probabilities::Regions2d ) );
3531 break;
3532 case 3 :
3533 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) );
3534 break;
3535 case 4 :
3537 break;
3538 case 5 :
3540 break;
3541 case 6 :
3542 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) );
3543 break;
3544 case 7 :
3546 break;
3547 case 8 :
3549 break;
3550 case 9 :
3552 break;
3553 case 10 :
3555 break;
3556 case 11 :
3557 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) );
3558 break;
3559 case 12 :
3561 break;
3562 }
3563 }
3564
3565 if( a_probability2d != nullptr ) {
3566 switch( type ) {
3567 case 0 :
3568 break;
3569 case 1 :
3570 static_cast<Probabilities::XYs2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3571 break;
3572 case 2 :
3573 static_cast<Probabilities::Regions2d * >( a_probability2d )->serialize( a_buffer, a_mode );
3574 break;
3575 case 3 :
3576 static_cast<Probabilities::Isotropic2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3577 break;
3578 case 4 :
3579 static_cast<Probabilities::DiscreteGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3580 break;
3581 case 5 :
3582 static_cast<Probabilities::PrimaryGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3583 break;
3584 case 6 :
3585 static_cast<Probabilities::Recoil2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3586 break;
3587 case 7 :
3588 static_cast<Probabilities::NBodyPhaseSpace2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3589 break;
3590 case 8 :
3591 static_cast<Probabilities::Evaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3592 break;
3593 case 9 :
3594 static_cast<Probabilities::GeneralEvaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3595 break;
3596 case 10 :
3597 static_cast<Probabilities::SimpleMaxwellianFission2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3598 break;
3599 case 11 :
3600 static_cast<Probabilities::Watt2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3601 break;
3602 case 12 :
3603 static_cast<Probabilities::WeightedFunctionals2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3604 break;
3605 }
3606 }
3607
3608 return( a_probability2d );
3609}
LUPI_HOST_DEVICE ProbabilityBase2dType ProbabilityBase2dClass(Probabilities::ProbabilityBase2d *funct)

Referenced by MCGIDI::Distributions::Uncorrelated::serialize(), and MCGIDI::HeatedReactionCrossSectionContinuousEnergy::serialize().

◆ serializeProbability2d_d1()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * MCGIDI::serializeProbability2d_d1 ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Probabilities::ProbabilityBase2d_d1 * a_probability2d )

Definition at line 3614 of file MCGIDI_functions.cc.

3615 {
3616
3617 int type = 0;
3618
3619 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3620 ProbabilityBase2dType pType = ProbabilityBase2dClass( a_probability2d );
3621
3622 switch( pType ) {
3624 break;
3626 type = 1;
3627 break;
3629 type = 2;
3630 break;
3632 type = 3;
3633 break;
3635 type = 4;
3636 break;
3638 type = 5;
3639 break;
3641 type = 6;
3642 break;
3644 type = 7;
3645 break;
3647 type = 8;
3648 break;
3650 type = 9;
3651 break;
3653 type = 10;
3654 break;
3656 type = 11;
3657 break;
3658 default:
3659 LUPI_THROW( "Probabilities::ProbabilityBase2d_d1: Unsupported ProbabilityBase2d_d1." );
3660 }
3661 }
3662
3663 DATA_MEMBER_INT( type, a_buffer, a_mode );
3664
3665 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3666 switch( type ) {
3667 case 0 :
3668 a_probability2d = nullptr;
3669 break;
3670 case 1 :
3671 if( a_buffer.m_placement != nullptr ) {
3672 a_probability2d = new(a_buffer.m_placement) Probabilities::XYs2d;
3673 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) ); }
3674 else {
3675 a_probability2d = new Probabilities::XYs2d;
3676 }
3677 break;
3678 case 2 :
3679 if( a_buffer.m_placement != nullptr ) {
3680 a_probability2d = new(a_buffer.m_placement) Probabilities::Regions2d;
3681 a_buffer.incrementPlacement( sizeof( Probabilities::Regions2d ) ); }
3682 else {
3683 a_probability2d = new Probabilities::Regions2d;
3684 }
3685 break;
3686 case 3 :
3687 if( a_buffer.m_placement != nullptr ) {
3688 a_probability2d = new(a_buffer.m_placement) Probabilities::Isotropic2d;
3689 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) ); }
3690 else {
3691 a_probability2d = new Probabilities::Isotropic2d;
3692 }
3693 break;
3694 case 4 :
3695 if( a_buffer.m_placement != nullptr ) {
3696 a_probability2d = new(a_buffer.m_placement) Probabilities::DiscreteGamma2d;
3697 a_buffer.incrementPlacement( sizeof( Probabilities::DiscreteGamma2d ) ); }
3698 else {
3699 a_probability2d = new Probabilities::DiscreteGamma2d;
3700 }
3701 break;
3702 case 5 :
3703 if( a_buffer.m_placement != nullptr ) {
3704 a_probability2d = new(a_buffer.m_placement) Probabilities::PrimaryGamma2d;
3705 a_buffer.incrementPlacement( sizeof( Probabilities::PrimaryGamma2d ) ); }
3706 else {
3707 a_probability2d = new Probabilities::PrimaryGamma2d;
3708 }
3709 break;
3710 case 6 :
3711 if( a_buffer.m_placement != nullptr ) {
3712 a_probability2d = new(a_buffer.m_placement) Probabilities::Recoil2d;
3713 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) ); }
3714 else {
3715 a_probability2d = new Probabilities::Recoil2d;
3716 }
3717 break;
3718 case 7 :
3719 if( a_buffer.m_placement != nullptr ) {
3720 a_probability2d = new(a_buffer.m_placement) Probabilities::NBodyPhaseSpace2d;
3722 else {
3723 a_probability2d = new Probabilities::NBodyPhaseSpace2d;
3724 }
3725 break;
3726 case 8 :
3727 if( a_buffer.m_placement != nullptr ) {
3728 a_probability2d = new(a_buffer.m_placement) Probabilities::Evaporation2d;
3729 a_buffer.incrementPlacement( sizeof( Probabilities::Evaporation2d ) ); }
3730 else {
3731 a_probability2d = new Probabilities::Evaporation2d;
3732 }
3733 break;
3734 case 9 :
3735 if( a_buffer.m_placement != nullptr ) {
3736 a_probability2d = new(a_buffer.m_placement) Probabilities::GeneralEvaporation2d;
3738 else {
3739 a_probability2d = new Probabilities::GeneralEvaporation2d;
3740 }
3741 break;
3742 case 10 :
3743 if( a_buffer.m_placement != nullptr ) {
3744 a_probability2d = new(a_buffer.m_placement) Probabilities::SimpleMaxwellianFission2d;
3746 else {
3747 a_probability2d = new Probabilities::SimpleMaxwellianFission2d;
3748 }
3749 break;
3750 case 11 :
3751 if( a_buffer.m_placement != nullptr ) {
3752 a_probability2d = new(a_buffer.m_placement) Probabilities::Watt2d;
3753 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) ); }
3754 else {
3755 a_probability2d = new Probabilities::Watt2d;
3756 }
3757 break;
3758 }
3759 }
3760
3761 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3762 switch( type ) {
3763 case 0 :
3764 break;
3765 case 1 :
3766 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) );
3767 break;
3768 case 2 :
3769 a_buffer.incrementPlacement( sizeof( Probabilities::Regions2d ) );
3770 break;
3771 case 3 :
3772 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) );
3773 break;
3774 case 4 :
3776 break;
3777 case 5 :
3779 break;
3780 case 6 :
3781 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) );
3782 break;
3783 case 7 :
3785 break;
3786 case 8 :
3788 break;
3789 case 9 :
3791 break;
3792 case 10 :
3794 break;
3795 case 11 :
3796 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) );
3797 break;
3798 }
3799 }
3800
3801 if( a_probability2d != nullptr ) {
3802 switch( type ) {
3803 case 0 :
3804 break;
3805 case 1 :
3806 static_cast<Probabilities::XYs2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3807 break;
3808 case 2 :
3809 static_cast<Probabilities::Regions2d * >( a_probability2d )->serialize( a_buffer, a_mode );
3810 break;
3811 case 3 :
3812 static_cast<Probabilities::Isotropic2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3813 break;
3814 case 4 :
3815 static_cast<Probabilities::DiscreteGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3816 break;
3817 case 5 :
3818 static_cast<Probabilities::PrimaryGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3819 break;
3820 case 6 :
3821 static_cast<Probabilities::Recoil2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3822 break;
3823 case 7 :
3824 static_cast<Probabilities::NBodyPhaseSpace2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3825 break;
3826 case 8 :
3827 static_cast<Probabilities::Evaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3828 break;
3829 case 9 :
3830 static_cast<Probabilities::GeneralEvaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3831 break;
3832 case 10 :
3833 static_cast<Probabilities::SimpleMaxwellianFission2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3834 break;
3835 case 11 :
3836 static_cast<Probabilities::Watt2d *>( a_probability2d )->serialize( a_buffer, a_mode );
3837 break;
3838 }
3839 }
3840
3841 return( a_probability2d );
3842}

Referenced by MCGIDI::Distributions::AngularEnergyMC::serialize(), MCGIDI::Distributions::AngularTwoBody::serialize(), MCGIDI::Distributions::EnergyAngularMC::serialize(), MCGIDI::Distributions::KalbachMann::serialize(), MCGIDI::Distributions::Uncorrelated::serialize(), MCGIDI::Probabilities::WeightedFunctionals2d::serialize(), and MCGIDI::Probabilities::XYs3d::serialize().

◆ serializeProbability2d_d2()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d2 * MCGIDI::serializeProbability2d_d2 ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Probabilities::ProbabilityBase2d_d2 * a_probability2d )

Definition at line 3847 of file MCGIDI_functions.cc.

3848 {
3849
3850 int type = 0;
3851
3852 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
3853 ProbabilityBase2dType pType = ProbabilityBase2dClass( a_probability2d );
3854
3855 switch( pType ) {
3857 break;
3859 type = 1;
3860 break;
3862 type = 3;
3863 break;
3865 type = 4;
3866 break;
3868 type = 5;
3869 break;
3871 type = 6;
3872 break;
3874 type = 7;
3875 break;
3877 type = 8;
3878 break;
3880 type = 9;
3881 break;
3883 type = 10;
3884 break;
3886 type = 11;
3887 break;
3888 default:
3889 LUPI_THROW( "Probabilities::ProbabilityBase2d_d1: Unsupported ProbabilityBase2d_d2" );
3890 }
3891 }
3892
3893 DATA_MEMBER_INT( type, a_buffer, a_mode );
3894
3895 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
3896 switch( type ) {
3897 case 0 :
3898 a_probability2d = nullptr;
3899 break;
3900 case 1 :
3901 if( a_buffer.m_placement != nullptr ) {
3902 a_probability2d = new(a_buffer.m_placement) Probabilities::XYs2d;
3903 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) ); }
3904 else {
3905 a_probability2d = new Probabilities::XYs2d;
3906 }
3907 break;
3908 case 3 :
3909 if( a_buffer.m_placement != nullptr ) {
3910 a_probability2d = new(a_buffer.m_placement) Probabilities::Isotropic2d;
3911 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) ); }
3912 else {
3913 a_probability2d = new Probabilities::Isotropic2d;
3914 }
3915 break;
3916 case 4 :
3917 if( a_buffer.m_placement != nullptr ) {
3918 a_probability2d = new(a_buffer.m_placement) Probabilities::DiscreteGamma2d;
3919 a_buffer.incrementPlacement( sizeof( Probabilities::DiscreteGamma2d ) ); }
3920 else {
3921 a_probability2d = new Probabilities::DiscreteGamma2d;
3922 }
3923 break;
3924 case 5 :
3925 if( a_buffer.m_placement != nullptr ) {
3926 a_probability2d = new(a_buffer.m_placement) Probabilities::PrimaryGamma2d;
3927 a_buffer.incrementPlacement( sizeof( Probabilities::PrimaryGamma2d ) ); }
3928 else {
3929 a_probability2d = new Probabilities::PrimaryGamma2d;
3930 }
3931 break;
3932 case 6 :
3933 if( a_buffer.m_placement != nullptr ) {
3934 a_probability2d = new(a_buffer.m_placement) Probabilities::Recoil2d;
3935 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) ); }
3936 else {
3937 a_probability2d = new Probabilities::Recoil2d;
3938 }
3939 break;
3940 case 7 :
3941 if( a_buffer.m_placement != nullptr ) {
3942 a_probability2d = new(a_buffer.m_placement) Probabilities::NBodyPhaseSpace2d;
3944 else {
3945 a_probability2d = new Probabilities::NBodyPhaseSpace2d;
3946 }
3947 break;
3948 case 8 :
3949 if( a_buffer.m_placement != nullptr ) {
3950 a_probability2d = new(a_buffer.m_placement) Probabilities::Evaporation2d;
3951 a_buffer.incrementPlacement( sizeof( Probabilities::Evaporation2d ) ); }
3952 else {
3953 a_probability2d = new Probabilities::Evaporation2d;
3954 }
3955 break;
3956 case 9 :
3957 if( a_buffer.m_placement != nullptr ) {
3958 a_probability2d = new(a_buffer.m_placement) Probabilities::GeneralEvaporation2d;
3960 else {
3961 a_probability2d = new Probabilities::GeneralEvaporation2d;
3962 }
3963 break;
3964 case 10 :
3965 if( a_buffer.m_placement != nullptr ) {
3966 a_probability2d = new(a_buffer.m_placement) Probabilities::SimpleMaxwellianFission2d;
3968 else {
3969 a_probability2d = new Probabilities::SimpleMaxwellianFission2d;
3970 }
3971 break;
3972 case 11 :
3973 if( a_buffer.m_placement != nullptr ) {
3974 a_probability2d = new(a_buffer.m_placement) Probabilities::Watt2d;
3975 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) ); }
3976 else {
3977 a_probability2d = new Probabilities::Watt2d;
3978 }
3979 break;
3980 }
3981 }
3982
3983 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
3984 switch( type ) {
3985 case 0 :
3986 break;
3987 case 1 :
3988 a_buffer.incrementPlacement( sizeof( Probabilities::XYs2d ) );
3989 break;
3990 case 3 :
3991 a_buffer.incrementPlacement( sizeof( Probabilities::Isotropic2d ) );
3992 break;
3993 case 4 :
3995 break;
3996 case 5 :
3998 break;
3999 case 6 :
4000 a_buffer.incrementPlacement( sizeof( Probabilities::Recoil2d ) );
4001 break;
4002 case 7 :
4004 break;
4005 case 8 :
4007 break;
4008 case 9 :
4010 break;
4011 case 10 :
4013 break;
4014 case 11 :
4015 a_buffer.incrementPlacement( sizeof( Probabilities::Watt2d ) );
4016 break;
4017 }
4018 }
4019
4020 if( a_probability2d != nullptr ) {
4021 switch( type ) {
4022 case 0 :
4023 break;
4024 case 1 :
4025 static_cast<Probabilities::XYs2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4026 break;
4027 case 3 :
4028 static_cast<Probabilities::Isotropic2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4029 break;
4030 case 4 :
4031 static_cast<Probabilities::DiscreteGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4032 break;
4033 case 5 :
4034 static_cast<Probabilities::PrimaryGamma2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4035 break;
4036 case 6 :
4037 static_cast<Probabilities::Recoil2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4038 break;
4039 case 7 :
4040 static_cast<Probabilities::NBodyPhaseSpace2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4041 break;
4042 case 8 :
4043 static_cast<Probabilities::Evaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4044 break;
4045 case 9 :
4046 static_cast<Probabilities::GeneralEvaporation2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4047 break;
4048 case 10 :
4049 static_cast<Probabilities::SimpleMaxwellianFission2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4050 break;
4051 case 11 :
4052 static_cast<Probabilities::Watt2d *>( a_probability2d )->serialize( a_buffer, a_mode );
4053 break;
4054 }
4055 }
4056
4057 return( a_probability2d );
4058}

Referenced by MCGIDI::Probabilities::Regions2d::serialize().

◆ serializeProbability3d()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * MCGIDI::serializeProbability3d ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Probabilities::ProbabilityBase3d * a_probability3d )

Definition at line 4072 of file MCGIDI_functions.cc.

4072 {
4073
4074 int type = 0;
4075
4076 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
4077 ProbabilityBase3dType pType = ProbabilityBase3dClass( a_probability3d );
4078
4079 switch( pType ) {
4081 break;
4083 type = 1;
4084 break;
4085 }
4086 }
4087
4088 DATA_MEMBER_INT( type, a_buffer, a_mode );
4089
4090 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
4091 switch( type ) {
4092 case 0 :
4093 a_probability3d = nullptr;
4094 break;
4095 case 1 :
4096 if( a_buffer.m_placement != nullptr ) {
4097 a_probability3d = new(a_buffer.m_placement) Probabilities::XYs3d;
4098 a_buffer.incrementPlacement( sizeof( Probabilities::XYs3d ) ); }
4099 else {
4100 a_probability3d = new Probabilities::XYs3d;
4101 }
4102 break;
4103 }
4104 }
4105 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
4106 switch( type ) {
4107 case 0 :
4108 break;
4109 case 1 :
4110 a_buffer.incrementPlacement( sizeof( Probabilities::XYs3d ) );
4111 break;
4112 }
4113 }
4114
4115 if( a_probability3d != nullptr ) {
4116 switch( type ) {
4117 case 0 :
4118 break;
4119 case 1 :
4120 static_cast<Probabilities::XYs3d *>( a_probability3d )->serialize( a_buffer, a_mode );
4121 break;
4122 }
4123 }
4124
4125 return( a_probability3d );
4126}
LUPI_HOST_DEVICE ProbabilityBase3dType ProbabilityBase3dClass(Probabilities::ProbabilityBase3d *funct)

Referenced by MCGIDI::Distributions::AngularEnergyMC::serialize(), and MCGIDI::Distributions::EnergyAngularMC::serialize().

◆ serializeProducts()

LUPI_HOST_DEVICE void MCGIDI::serializeProducts ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Vector< Product * > & a_products )

This method serializes a_products for broadcasting as needed for MPI and GPUs. The method can count the number of required bytes, pack a_products or unpack a_products 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.
a_products[in] The products to serialize.

Definition at line 367 of file MCGIDI_misc.cc.

367 {
368
369 std::size_t vectorSize = a_products.size( );
370 int vectorSizeInt = (int) vectorSize;
371 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
372 vectorSize = (std::size_t) vectorSizeInt;
373
374 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
375 a_products.resize( vectorSize, &a_buffer.m_placement );
376 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
377 if( a_buffer.m_placement != nullptr ) {
378 a_products[vectorIndex] = new(a_buffer.m_placement) Product;
379 a_buffer.incrementPlacement( sizeof( Product ) );
380 }
381 else {
382 a_products[vectorIndex] = new Product;
383 }
384 } }
385 else if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
386 a_buffer.m_placement += a_products.internalSize( );
387 a_buffer.incrementPlacement( sizeof( Product ) * vectorSize );
388 }
389
390 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
391 a_products[vectorIndex]->serialize( a_buffer, a_mode );
392 }
393}

Referenced by MCGIDI::OutputChannel::serialize(), and MCGIDI::Reaction::serialize().

◆ serializeQs()

LUPI_HOST_DEVICE void MCGIDI::serializeQs ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode,
Vector< Functions::Function1d_d1 * > & a_Qs )

This method serializes a_Qs for broadcasting as needed for MPI and GPUs. The method can count the number of required bytes, pack a_Qs or unpack a_Qs 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.
a_Qs[in] The Q functions to serialize.

Definition at line 441 of file MCGIDI_misc.cc.

441 {
442
443 std::size_t vectorSize = a_Qs.size( );
444 int vectorSizeInt = (int) vectorSize;
445 DATA_MEMBER_INT( vectorSizeInt, a_buffer, a_mode );
446 vectorSize = (std::size_t) vectorSizeInt;
447
448 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
449 a_Qs.resize( vectorSize, &a_buffer.m_placement ); }
450 else if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
451 a_buffer.m_placement += a_Qs.internalSize( );
452 }
453
454 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
455 a_Qs[vectorIndex] = serializeFunction1d_d1( a_buffer, a_mode, a_Qs[vectorIndex] );
456 }
457}
LUPI_HOST_DEVICE Functions::Function1d_d1 * serializeFunction1d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d)

Referenced by MCGIDI::Reaction::serialize().

◆ serializeURR_mode()

LUPI_HOST_DEVICE Transporting::URR_mode MCGIDI::serializeURR_mode ( Transporting::URR_mode a_URR_mode,
LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )

This method serializes a Transporting::URR_mode value.

Parameters
a_URR_mode[in] The inputted Transporting::URR_mode value.
a_buffer[in] The buffer to read or write data to depending on a_mode.
a_mode[in] Specifies the action of this method.
Returns
The Transporting::URR_mode value.

Definition at line 349 of file MCGIDI_URR.cc.

349 {
350
351 int type = 0;
352 switch( a_URR_mode ) {
354 break;
356 type = 1;
357 break;
359 type = 2;
360 break;
361 }
362 DATA_MEMBER_INT( type, a_buffer, a_mode );
363
364 if( type == 0 ) return( Transporting::URR_mode::none );
365 if( type == 1 ) return( Transporting::URR_mode::pdfs );
367}

Referenced by MCGIDI::HeatedCrossSectionContinuousEnergy::serialize(), and MCGIDI::HeatedReactionCrossSectionContinuousEnergy::serialize().

◆ upScatterModelABoostParticle()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::upScatterModelABoostParticle ( Sampling::Input & a_input,
RNG && a_rng,
Sampling::Product & a_product )
inline

This function boost a particle from one frame to another frame. The frames have a relative speed a_boostSpeed and cosine of angle a_boostMu between their z-axes. BRB FIXME, currently it is the x-axis.

Parameters
a_input[in] Instance containing a random number generator that returns a double in the range [0, 1).
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_product[in] The particle to boost.

Definition at line 2358 of file MCGIDI_headerSource.hpp.

2358 {
2359
2360 double C_rel = 1.0;
2361 if( a_input.m_relativeBeta != 0.0 ) {
2362 C_rel = ( a_input.m_projectileBeta - a_input.m_muLab * a_input.m_targetBeta ) / a_input.m_relativeBeta;
2363 if( C_rel > 1.0 ) C_rel = 1.0; // Handle round-off issue. Probably should check how big the issue is.
2364 if( C_rel < -1.0 ) C_rel = -1.0; // Handle round-off issue. Probably should check how big the issue is.
2365 }
2366 double S_rel = sqrt( 1.0 - C_rel * C_rel );
2367
2368 double pz_vz = a_product.m_pz_vz;
2369 a_product.m_pz_vz = C_rel * a_product.m_pz_vz + S_rel * a_product.m_px_vx;
2370 a_product.m_px_vx = -S_rel * pz_vz + C_rel * a_product.m_px_vx;
2371
2372 double targetSpeed = MCGIDI_speedOfLight_cm_sec * a_input.m_targetBeta;
2373 a_product.m_pz_vz += a_input.m_muLab * targetSpeed;
2374 a_product.m_px_vx += sqrt( 1.0 - a_input.m_muLab * a_input.m_muLab ) * targetSpeed;
2375
2376 double phi = 2.0 * M_PI * a_rng( );
2377 double sine = sin( phi );
2378 double cosine = cos( phi );
2379 double px_vx = a_product.m_px_vx;
2380 a_product.m_px_vx = cosine * a_product.m_px_vx - sine * a_product.m_py_vy;
2381 a_product.m_py_vy = sine * px_vx + cosine * a_product.m_py_vy;
2382
2383 double speed2 = a_product.m_px_vx * a_product.m_px_vx + a_product.m_py_vy * a_product.m_py_vy + a_product.m_pz_vz * a_product.m_pz_vz;
2385
2386 a_product.m_kineticEnergy = particleKineticEnergyFromBeta2( a_product.m_productMass, speed2 );
2387}
#define MCGIDI_speedOfLight_cm_sec
Definition MCGIDI.hpp:68
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE double particleKineticEnergyFromBeta2(double a_mass_unitOfEnergy, double a_particleBeta2)

Referenced by MCGIDI::Sampling::ProductHandler::add().

◆ vectorToSTD_vector() [1/2]

LUPI_HOST std::vector< double > MCGIDI::vectorToSTD_vector ( Vector< double > a_input)

This function returns a std::vector<double> that represents a_input.

Parameters
a_input[in] The input for the returned std::vector<double> instance.
Returns
A std::vector<double> instance.

Definition at line 510 of file MCGIDI_misc.cc.

510 {
511
512 std::vector<double> vector( a_input.size( ) );
513
514 std::size_t index = 0;
515 for( auto iter = a_input.begin( ); iter != a_input.end( ); ++iter, ++index ) vector[index] = *iter;
516
517 return( vector );
518}

Referenced by MCGIDI::HeatedCrossSectionContinuousEnergy::crossSectionAsGIDI_XYs1d(), and MCGIDI::HeatedReactionCrossSectionContinuousEnergy::crossSectionAsGIDI_XYs1d().

◆ vectorToSTD_vector() [2/2]

LUPI_HOST std::vector< double > MCGIDI::vectorToSTD_vector ( Vector< float > a_input)

This function returns a std::vector<double> that represents a_input.

Parameters
a_input[in] The input for the returned std::vector<double> instance.
Returns
A std::vector<double> instance.

Definition at line 528 of file MCGIDI_misc.cc.

528 {
529
530 std::vector<double> vector( a_input.size( ) );
531
532 std::size_t index = 0;
533 for( auto iter = a_input.begin( ); iter != a_input.end( ); ++iter, ++index ) vector[index] = *iter;
534
535 return( vector );
536}