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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::Uncorrelated:

Public Member Functions

LUPI_HOST_DEVICE Uncorrelated ()
LUPI_HOST Uncorrelated (GIDI::Distributions::Uncorrelated const &a_uncorrelated, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~Uncorrelated ()
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1angular () const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2denergy () 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 product for which the distribution is the product of uncorrelated angular (i.e., P(mu|E)) and energy (i.e., P(E'|E)) distributions.

Definition at line 102 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ Uncorrelated() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::Uncorrelated::Uncorrelated ( )

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

Definition at line 199 of file MCGIDI_distributions.cc.

199 :
200 m_angular( nullptr ),
201 m_energy( nullptr ) {
202
203}

◆ Uncorrelated() [2/2]

LUPI_HOST MCGIDI::Distributions::Uncorrelated::Uncorrelated ( GIDI::Distributions::Uncorrelated const & a_uncorrelated,
SetupInfo & a_setupInfo )
Parameters
a_uncorrelated[in] The GIDI::Distributions::Uncorrelated 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 210 of file MCGIDI_distributions.cc.

210 :
211 Distribution( Type::uncorrelated, a_uncorrelated, a_setupInfo ),
212 m_angular( Probabilities::parseProbability2d_d1( a_uncorrelated.angular( ), nullptr ) ),
213 m_energy( Probabilities::parseProbability2d( a_uncorrelated.energy( ), &a_setupInfo ) ) {
214
215}
LUPI_HOST ProbabilityBase2d_d1 * parseProbability2d_d1(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)
LUPI_HOST ProbabilityBase2d * parseProbability2d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite, SetupInfo *a_setupInfo)

◆ ~Uncorrelated()

LUPI_HOST_DEVICE MCGIDI::Distributions::Uncorrelated::~Uncorrelated ( )

Definition at line 220 of file MCGIDI_distributions.cc.

220 {
221
222 delete m_angular;
223 delete m_energy;
224}

Member Function Documentation

◆ angleBiasing() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::Uncorrelated::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 654 of file MCGIDI_headerSource.hpp.

655 {
656
657 if( productFrame( ) != GIDI::Frame::lab ) {
658 a_energy_out = 0.0;
659
660 double initialMass = projectileMass( ) + targetMass( );
661 double boostBeta = sqrt( a_energy_in * ( a_energy_in + 2. * projectileMass( ) ) ) / ( a_energy_in + initialMass ); // Good, even for projectileMass = 0.
662 double energy_out_com = m_energy->sample( a_energy_in, a_rng( ), a_rng );
663
664 if( productMass( ) == 0.0 ) {
665 double one_mu_beta = 1.0 - a_mu_lab * boostBeta;
666 double mu_com = ( a_mu_lab - boostBeta ) / one_mu_beta;
667 double Jacobian = ( 1.0 - boostBeta * boostBeta ) / ( one_mu_beta * one_mu_beta );
668
669 a_energy_out = sqrt( 1.0 - boostBeta * boostBeta ) * energy_out_com * ( 1.0 + mu_com * boostBeta );
670
671 return( Jacobian * m_angular->evaluate( a_energy_in, mu_com ) );
672 }
673
674 double productBeta = MCGIDI_particleBeta( productMass( ), energy_out_com );
675 double muPlus = 0.0, JacobianPlus = 0.0, muMinus = 0.0, JacobianMinus = 0.0;
676 int numberOfMus = muCOM_From_muLab( a_mu_lab, boostBeta, productBeta, muPlus, JacobianPlus, muMinus, JacobianMinus );
677
678 if( numberOfMus == 0 ) return( 0.0 );
679
680 double probability = JacobianPlus * m_angular->evaluate( a_energy_in, muPlus );
681
682 if( numberOfMus == 2 ) {
683 double probabilityMinus = JacobianMinus * m_angular->evaluate( a_energy_in, muMinus );
684
685 probability += probabilityMinus;
686 if( probabilityMinus > a_rng( ) * probability ) muPlus = muMinus;
687 }
688
689 double productBeta2 = productBeta * productBeta;
690 double productBetaLab2 = productBeta2 + boostBeta * boostBeta * ( 1.0 - productBeta2 * ( 1.0 - muPlus * muPlus ) ) + 2.0 * muPlus * productBeta * boostBeta;
691 productBetaLab2 /= 1.0 - muPlus * productBeta * boostBeta;
692 a_energy_out = MCGIDI::particleKineticEnergyFromBeta2( productMass( ), productBetaLab2 );
693
694 return( probability );
695 }
696
697 a_energy_out = m_energy->sample( a_energy_in, a_rng( ), a_rng );
698 return( m_angular->evaluate( a_energy_in, a_mu_lab ) );
699}
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71
LUPI_HOST_DEVICE double targetMass() const
LUPI_HOST_DEVICE GIDI::Frame productFrame() const
LUPI_HOST_DEVICE double productMass() const
LUPI_HOST_DEVICE double projectileMass() const
LUPI_HOST_DEVICE double particleKineticEnergyFromBeta2(double a_mass_unitOfEnergy, double a_particleBeta2)
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)

◆ angleBiasing() [2/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::Uncorrelated::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::Uncorrelated::angular ( ) const
inline

Returns the value of the m_angular.

Definition at line 113 of file MCGIDI_distributions.hpp.

Referenced by Uncorrelated().

◆ energy()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d * MCGIDI::Distributions::Uncorrelated::energy ( ) const
inline

Returns the value of the m_energy.

Definition at line 114 of file MCGIDI_distributions.hpp.

Referenced by Uncorrelated().

◆ sample()

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

This method samples the outgoing product data by sampling the outgoing energy E' and mu from the uncorrelated E and mu probabilities. 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 630 of file MCGIDI_headerSource.hpp.

630 {
631
632 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
633 a_input.m_mu = m_angular->sample( a_X, a_rng( ), a_rng );
634 a_input.m_energyOut1 = m_energy->sample( a_X, a_rng( ), a_rng );
635 a_input.m_phi = 2. * M_PI * a_rng( );
636 a_input.m_frame = productFrame( );
637}
#define M_PI
Definition SbMath.h:33

◆ serialize()

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

234 {
235
236 Distribution::serialize( a_buffer, a_mode );
237 m_angular = serializeProbability2d_d1( a_buffer, a_mode, m_angular );
238 m_energy = serializeProbability2d( a_buffer, a_mode, m_energy );
239}
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d * serializeProbability2d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d *a_probability2d)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * serializeProbability2d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d1 *a_probability2d)

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