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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::AngularEnergyMC:

Public Member Functions

LUPI_HOST_DEVICE AngularEnergyMC ()
LUPI_HOST AngularEnergyMC (GIDI::Distributions::AngularEnergyMC const &a_angularEnergyMC, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~AngularEnergyMC ()
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1angular () const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3denergyGivenAngular () const
template<typename RNG>
LUPI_HOST_DEVICE void sample (double a_X, 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, double a_mu_lab, 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 particle where the distribution is give as P(mu|E) * P(E'|E,mu) where E is the projectile's energy, E' is the product's outgoing energy, mu is the cosine of the product's outgoing angle relative to the projectile's velocity, P(mu|E) is the probability for mu given E and (P(E'|E,mu) is the probability for E' given E and mu.

Definition at line 177 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ AngularEnergyMC() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::AngularEnergyMC::AngularEnergyMC ( )

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

Definition at line 360 of file MCGIDI_distributions.cc.

360 :
361 m_angular( nullptr ),
362 m_energyGivenAngular( nullptr ) {
363
364}

◆ AngularEnergyMC() [2/2]

LUPI_HOST MCGIDI::Distributions::AngularEnergyMC::AngularEnergyMC ( GIDI::Distributions::AngularEnergyMC const & a_angularEnergyMC,
SetupInfo & a_setupInfo )
Parameters
a_angularEnergyMC[in] The GIDI::Distributions::AngularEnergyMC instance whose data is to be used to construct this.
a_setupInfo[in] Used internally when constructing a Protare to pass information to other constructors.

Definition at line 371 of file MCGIDI_distributions.cc.

371 :
372 Distribution( Type::angularEnergyMC, a_angularEnergyMC, a_setupInfo ),
373 m_angular( Probabilities::parseProbability2d_d1( a_angularEnergyMC.angular( ), nullptr ) ),
374 m_energyGivenAngular( Probabilities::parseProbability3d( a_angularEnergyMC.angularEnergy( ) ) ) {
375
376}
LUPI_HOST ProbabilityBase3d * parseProbability3d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST ProbabilityBase2d_d1 * parseProbability2d_d1(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)

◆ ~AngularEnergyMC()

LUPI_HOST_DEVICE MCGIDI::Distributions::AngularEnergyMC::~AngularEnergyMC ( )

Definition at line 381 of file MCGIDI_distributions.cc.

381 {
382
383 delete m_angular;
384 delete m_energyGivenAngular;
385}

Member Function Documentation

◆ angleBiasing() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::AngularEnergyMC::angleBiasing ( LUPI_maybeUnused Reaction const * a_reaction,
LUPI_maybeUnused double a_temperature,
double a_energy_in,
double a_mu_lab,
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 858 of file MCGIDI_headerSource.hpp.

859 {
860
861 if( productFrame( ) != GIDI::Frame::lab ) LUPI_THROW( "AngularEnergyMC::angleBiasing: center-of-mass not supported." );
862
863 a_energy_out = m_energyGivenAngular->sample( a_energy_in, a_mu_lab, a_mu_lab, a_rng( ), a_rng );
864 return( m_angular->evaluate( a_energy_in, a_mu_lab ) );
865}
#define LUPI_THROW(arg)
LUPI_HOST_DEVICE GIDI::Frame productFrame() const

◆ angleBiasing() [2/2]

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

◆ angular()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * MCGIDI::Distributions::AngularEnergyMC::angular ( ) const
inline

Returns the value of the m_angular.

Definition at line 188 of file MCGIDI_distributions.hpp.

Referenced by AngularEnergyMC().

◆ energyGivenAngular()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * MCGIDI::Distributions::AngularEnergyMC::energyGivenAngular ( ) const
inline

Returns the value of the m_energyGivenAngular.

Definition at line 189 of file MCGIDI_distributions.hpp.

◆ sample()

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

This method samples the outgoing product data by sampling the outgoing mu from the probability P(mu|E) and then samples E' from the probability P(E'|E,mu). It also samples the outgoing phi uniformly between 0 and 2 pi.

Parameters
a_X[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 832 of file MCGIDI_headerSource.hpp.

832 {
833
834 double mu_1, mu_2;
835
836 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
837 a_input.m_mu = m_angular->sample2dOf3d( a_X, a_rng( ), a_rng, &mu_1, &mu_2 );
838 a_input.m_energyOut1 = m_energyGivenAngular->sample( a_X, mu_1, mu_2, a_rng( ), a_rng );
839 a_input.m_phi = 2. * M_PI * a_rng( );
840 a_input.m_frame = productFrame( );
841}
#define M_PI
Definition SbMath.h:33

◆ serialize()

LUPI_HOST_DEVICE void MCGIDI::Distributions::AngularEnergyMC::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 395 of file MCGIDI_distributions.cc.

395 {
396
397 Distribution::serialize( a_buffer, a_mode );
398 m_angular = serializeProbability2d_d1( a_buffer, a_mode, m_angular );
399 m_energyGivenAngular = serializeProbability3d( a_buffer, a_mode, m_energyGivenAngular );
400}
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * serializeProbability2d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d1 *a_probability2d)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * serializeProbability3d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase3d *a_probability3d)

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