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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::CoherentElasticTNSL:

Public Member Functions

LUPI_HOST_DEVICE CoherentElasticTNSL ()
LUPI_HOST CoherentElasticTNSL (GIDI::DoubleDifferentialCrossSection::n_ThermalNeutronScatteringLaw::CoherentElastic const *a_coherentElasticTNSL, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~CoherentElasticTNSL ()
template<typename RNG>
LUPI_HOST_DEVICE void sample (double a_energy, Sampling::Input &a_input, RNG &&a_rng) const
template<typename RNG>
LUPI_HOST_DEVICE double angleBiasing (Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
template<typename RNG>
LUPI_HOST_DEVICE double angleBiasing (LUPI_maybeUnused Reaction const *a_reaction, LUPI_maybeUnused double a_temperature, double a_energy_in, LUPI_maybeUnused double a_mu_lab, LUPI_maybeUnused RNG &&a_rng, double &a_energy_out) const
Public Member Functions inherited from MCGIDI::Distributions::Distribution
LUPI_HOST_DEVICE Distribution ()
LUPI_HOST Distribution (Type a_type, GIDI::Distributions::Distribution const &a_distribution, SetupInfo &a_setupInfo)
LUPI_HOST Distribution (Type a_type, GIDI::Frame a_productFrame, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION ~Distribution () MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE Type type () const
LUPI_HOST_DEVICE GIDI::Frame productFrame () const
LUPI_HOST_DEVICE double projectileMass () const
LUPI_HOST_DEVICE double targetMass () const
LUPI_HOST_DEVICE double productMass () const
LUPI_HOST void setModelDBRC_data (Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
template<typename RNG>
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION void sample (double a_X, Sampling::Input &a_input, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
template<typename RNG>
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double angleBiasing (Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
template<typename RNG>
LUPI_HOST_DEVICE void sample (double a_X, MCGIDI::Sampling::Input &a_input, RNG &&a_rng) const
template<typename RNG>
LUPI_HOST_DEVICE double angleBiasing (Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const

Detailed Description

This class represents the distribution for an outgoing product whose distribution is TNSL coherent elastic scattering. This class samples directly from the Debye/Waller function.

Definition at line 382 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ CoherentElasticTNSL() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::CoherentElasticTNSL::CoherentElasticTNSL ( )

Constructor for the CoherentElasticTNSL class.

Definition at line 1225 of file MCGIDI_distributions.cc.

1225 :
1226 m_temperatureInterpolation( Interpolation::LINLIN ) {
1227
1228}

◆ CoherentElasticTNSL() [2/2]

LUPI_HOST MCGIDI::Distributions::CoherentElasticTNSL::CoherentElasticTNSL ( GIDI::DoubleDifferentialCrossSection::n_ThermalNeutronScatteringLaw::CoherentElastic const * a_coherentElasticTNSL,
SetupInfo & a_setupInfo )

Constructor for the CoherentElasticTNSL class.

Parameters
a_coherentElasticTNSL[in] GIDI::CoherentElastic instance containing the Debye/Waller data.
a_setupInfo[in] Used internally when constructing a Protare to pass information to other constructors.

Definition at line 1237 of file MCGIDI_distributions.cc.

1238 :
1240 m_temperatureInterpolation( Interpolation::LINLIN ) {
1241
1242 GIDI::DoubleDifferentialCrossSection::n_ThermalNeutronScatteringLaw::S_table const &s_table = a_coherentElasticTNSL->s_table( );
1243 GIDI::Functions::Gridded2d const *gridded2d = dynamic_cast<GIDI::Functions::Gridded2d const *>( s_table.function2d( ) );
1244 GIDI::Axes const &axes = gridded2d->axes( );
1245
1246 GIDI::Grid const *axis = dynamic_cast<GIDI::Grid const *>( axes[0] );
1247 m_temperatureInterpolation = GIDI2MCGIDI_interpolation( ptwXY_stringToInterpolation( axis->interpolation( ).c_str( ) ) );
1248
1249 double temperatureToMeV_K = 1.0;
1250 if( axis->unit( ) == "K" ) temperatureToMeV_K = 8.617330337217212e-11; // This is a kludge until units are properly supported.
1251 nf_Buffer<double> const *grid = &axis->values( );
1252 m_temperatures.resize( grid->size( ) );
1253 for( std::size_t index = 0; index < grid->size( ); ++index ) m_temperatures[index] = temperatureToMeV_K * (*grid)[index];
1254
1255 axis = dynamic_cast<GIDI::Grid const *>( axes[1] );
1256 grid = &axis->values( );
1257 m_energies.resize( grid->size( ) );
1258 for( std::size_t index = 0; index < grid->size( ); ++index ) m_energies[index] = (*grid)[index];
1259
1260 GIDI::Array::FullArray fullArray = gridded2d->array( ).constructArray( );
1261 m_S_table.resize( fullArray.size( ) );
1262 for( std::size_t index = 0; index < fullArray.m_flattenedValues.size( ); ++index ) m_S_table[index] = fullArray.m_flattenedValues[index];
1263}
std::size_t size() const
Definition GIDI.hpp:900
std::vector< double > m_flattenedValues
Definition GIDI.hpp:898
LUPI_HOST_DEVICE Interpolation GIDI2MCGIDI_interpolation(ptwXY_interpolation a_interpolation)
ptwXY_interpolation ptwXY_stringToInterpolation(char const *interpolationString)
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668

◆ ~CoherentElasticTNSL()

LUPI_HOST_DEVICE MCGIDI::Distributions::CoherentElasticTNSL::~CoherentElasticTNSL ( )
inline

Definition at line 394 of file MCGIDI_distributions.hpp.

394{}

Member Function Documentation

◆ angleBiasing() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::CoherentElasticTNSL::angleBiasing ( LUPI_maybeUnused Reaction const * a_reaction,
LUPI_maybeUnused double a_temperature,
double a_energy_in,
LUPI_maybeUnused double a_mu_lab,
LUPI_maybeUnused RNG && a_rng,
double & a_energy_out ) const

Returns the probability for a projectile with energy a_energy_in to cause a particle to be emitted at angle a_mu_lab as seen in the lab frame. a_energy_out is the sampled outgoing energy.

Parameters
a_reaction[in] The reaction containing the particle which this distribution describes.
a_temperature[in] The temperature of the material.
a_energy_in[in] The energy of the incident particle.
a_mu_lab[in] The desired mu in the lab frame for the emitted particle.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_energy_out[in] The energy of the emitted outgoing particle.
Returns
The probability of emitting outgoing particle into lab angle a_mu_lab.

Definition at line 1475 of file MCGIDI_headerSource.hpp.

1476 {
1477
1478 double probability = 0.0;
1479
1480 a_energy_out = a_energy_in;
1481
1482 return( probability );
1483}

◆ angleBiasing() [2/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::CoherentElasticTNSL::angleBiasing ( Reaction const * a_reaction,
double a_temperature,
double a_energy_in,
double a_mu_lab,
RNG && a_rng,
double & a_energy_out ) const

◆ sample()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::Distributions::CoherentElasticTNSL::sample ( double a_energy,
Sampling::Input & a_input,
RNG && a_rng ) const

This method samples the outgoing neutron data for coherent elastic TSNL from the Debye/Waller function.

Parameters
a_energy[in] The energy of the projectile.
a_input[in] Sample options requested by user.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).

Definition at line 1423 of file MCGIDI_headerSource.hpp.

1424 {
1425
1426 if( a_energy <= m_energies[0] ) {
1427 a_input.m_mu = 1.0; }
1428 else {
1429 double temperature = a_input.temperature( );
1430 if( temperature < m_temperatures[0] ) temperature = m_temperatures[0];
1431 if( temperature > m_temperatures.back( ) ) temperature = m_temperatures.back( );
1432 std::size_t temperatureIndex = (std::size_t) MCGIDI::binarySearchVector( temperature, m_temperatures, true );
1433 double const *pointer1 = &m_S_table[temperatureIndex * m_energies.size( )];
1434
1435 double const *pointer2 = pointer1;
1436 double fractionFirstTemperature = 1.0;
1437 if( temperatureIndex != ( m_temperatures.size( ) - 1 ) ) {
1438 fractionFirstTemperature = ( m_temperatures[temperatureIndex+1] - temperature ) / ( m_temperatures[temperatureIndex+1] - m_temperatures[temperatureIndex] );
1439 pointer2 += m_energies.size( );
1440 }
1441 double fractionSecondTemperature = 1.0 - fractionFirstTemperature;
1442
1443 int intEnergyIndexMax = MCGIDI::binarySearchVector( a_energy, m_energies, true );
1444 std::size_t energyIndexMax = static_cast<std::size_t>( intEnergyIndexMax );
1445 if( a_energy == m_energies[energyIndexMax] ) --energyIndexMax;
1446
1447 double randomTotal = a_rng( ) * ( fractionFirstTemperature * pointer1[energyIndexMax] + fractionSecondTemperature * pointer2[energyIndexMax] );
1448 std::size_t energyIndex = 0;
1449 for( ; energyIndex < energyIndexMax; ++energyIndex ) {
1450 if( randomTotal <= fractionFirstTemperature * pointer1[energyIndex] + fractionSecondTemperature * pointer2[energyIndex] ) break;
1451 }
1452 a_input.m_mu = 1.0 - 2.0 * m_energies[energyIndex] / a_energy;
1453 }
1454
1455 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
1456 a_input.m_energyOut1 = a_energy;
1457 a_input.m_phi = 2.0 * M_PI * a_rng( );
1458}
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE int binarySearchVector(double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
Definition MCGIDI.hpp:318

◆ serialize()

LUPI_HOST_DEVICE void MCGIDI::Distributions::CoherentElasticTNSL::serialize ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )

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

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

Definition at line 1273 of file MCGIDI_distributions.cc.

1273 {
1274
1275 Distribution::serialize( a_buffer, a_mode );
1276
1277 DATA_MEMBER_VECTOR_DOUBLE( m_temperatures, a_buffer, a_mode );
1278 DATA_MEMBER_VECTOR_DOUBLE( m_energies, a_buffer, a_mode );
1279 DATA_MEMBER_VECTOR_DOUBLE( m_S_table, a_buffer, a_mode );
1280
1281 int interpolation = 0;
1282 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
1283 switch( m_temperatureInterpolation ) {
1284 case Interpolation::FLAT :
1285 break;
1287 interpolation = 1;
1288 break;
1290 interpolation = 2;
1291 break;
1293 interpolation = 3;
1294 break;
1296 interpolation = 4;
1297 break;
1299 interpolation = 5;
1300 break;
1301 }
1302 }
1303 DATA_MEMBER_INT( interpolation, a_buffer, a_mode );
1304 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
1305 switch( interpolation ) {
1306 case 0 :
1307 m_temperatureInterpolation = Interpolation::FLAT;
1308 break;
1309 case 1 :
1310 m_temperatureInterpolation = Interpolation::LINLIN;
1311 break;
1312 case 2 :
1313 m_temperatureInterpolation = Interpolation::LINLOG;
1314 break;
1315 case 3 :
1316 m_temperatureInterpolation = Interpolation::LOGLIN;
1317 break;
1318 case 4 :
1319 m_temperatureInterpolation = Interpolation::LOGLOG;
1320 break;
1321 case 5 :
1322 m_temperatureInterpolation = Interpolation::OTHER;
1323 break;
1324 }
1325 }
1326}
#define DATA_MEMBER_VECTOR_DOUBLE(member, buf, mode)
#define DATA_MEMBER_INT( member, buf, mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

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