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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::EnergyAngularMC:

Public Member Functions

LUPI_HOST_DEVICE EnergyAngularMC ()
LUPI_HOST EnergyAngularMC (GIDI::Distributions::EnergyAngularMC const &a_energyAngularMC, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~EnergyAngularMC ()
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1energy () const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3dangularGivenEnergy () 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(E'|E) * P(mu|E,E') 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(E'|E) is the probability for E' given E and (P(mu|E,E') is the probability for mu given E and E'.

Definition at line 151 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ EnergyAngularMC() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::EnergyAngularMC::EnergyAngularMC ( )

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

Definition at line 307 of file MCGIDI_distributions.cc.

307 :
308 m_energy( nullptr ),
309 m_angularGivenEnergy( nullptr ) {
310
311}

◆ EnergyAngularMC() [2/2]

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

318 :
319 Distribution( Type::energyAngularMC, a_energyAngularMC, a_setupInfo ),
320 m_energy( Probabilities::parseProbability2d_d1( a_energyAngularMC.energy( ), nullptr ) ),
321 m_angularGivenEnergy( Probabilities::parseProbability3d( a_energyAngularMC.energyAngular( ) ) ) {
322
323}
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)

◆ ~EnergyAngularMC()

LUPI_HOST_DEVICE MCGIDI::Distributions::EnergyAngularMC::~EnergyAngularMC ( )

Definition at line 328 of file MCGIDI_distributions.cc.

328 {
329
330 delete m_energy;
331 delete m_angularGivenEnergy;
332}

Member Function Documentation

◆ angleBiasing() [1/2]

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

774 {
775
776 double probability = 0.0;
777
779 a_energy_out = m_energy->sample( a_energy_in, a_rng( ), a_rng );
780
781 double initialMass = projectileMass( ) + targetMass( );
782 double boostBeta = sqrt( a_energy_in * ( a_energy_in + 2. * projectileMass( ) ) ) / ( a_energy_in + initialMass ); // Good, even for projectileMass = 0.
783 double energy_out_com = m_energy->sample( a_energy_in, a_rng( ), a_rng );
784
785 if( productMass( ) == 0.0 ) {
786 double one_mu_beta = 1.0 - a_mu_lab * boostBeta;
787 double mu_com = ( a_mu_lab - boostBeta ) / one_mu_beta;
788 double Jacobian = ( 1.0 - boostBeta * boostBeta ) / ( one_mu_beta * one_mu_beta );
789
790 a_energy_out = sqrt( 1.0 - boostBeta * boostBeta ) * energy_out_com * ( 1.0 + mu_com * boostBeta );
791
792 return( Jacobian * m_angularGivenEnergy->evaluate( a_energy_in, energy_out_com, mu_com ) );
793 }
794
795 double productBeta = MCGIDI_particleBeta( productMass( ), energy_out_com );
796 double muPlus = 0.0, JacobianPlus = 0.0, muMinus = 0.0, JacobianMinus = 0.0;
797 int numberOfMus = muCOM_From_muLab( a_mu_lab, boostBeta, productBeta, muPlus, JacobianPlus, muMinus, JacobianMinus );
798
799 if( numberOfMus == 0 ) return( 0.0 );
800
801 probability = JacobianPlus * m_angularGivenEnergy->evaluate( a_energy_in, energy_out_com, muPlus );
802
803 if( numberOfMus == 2 ) {
804 double probabilityMinus = JacobianMinus * m_angularGivenEnergy->evaluate( a_energy_in, energy_out_com, muMinus );
805
806 probability += probabilityMinus;
807 if( probabilityMinus > a_rng( ) * probability ) muPlus = muMinus;
808 }
809
810 double productBeta2 = productBeta * productBeta;
811 double productBetaLab2 = productBeta2 + boostBeta * boostBeta * ( 1.0 - productBeta2 * ( 1.0 - muPlus * muPlus ) ) + 2.0 * muPlus * productBeta * boostBeta;
812 productBetaLab2 /= 1.0 - muPlus * productBeta * boostBeta;
813 a_energy_out = MCGIDI::particleKineticEnergyFromBeta2( productMass( ), productBetaLab2 ); }
814 else {
815 a_energy_out = m_energy->sample( a_energy_in, a_rng( ), a_rng );
816 probability = m_angularGivenEnergy->evaluate( a_energy_in, a_energy_out, a_mu_lab );
817 }
818
819 return( probability );
820}
#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
@ centerOfMass
Definition GIDI.hpp:146
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::EnergyAngularMC::angleBiasing ( Reaction const * a_reaction,
double a_temperature,
double a_energy_in,
double a_mu_lab,
RNG && a_rng,
double & a_energy_out ) const

◆ angularGivenEnergy()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * MCGIDI::Distributions::EnergyAngularMC::angularGivenEnergy ( ) const
inline

Returns the value of the m_angularGivenEnergy.

Definition at line 163 of file MCGIDI_distributions.hpp.

◆ energy()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * MCGIDI::Distributions::EnergyAngularMC::energy ( ) const
inline

Returns the value of the m_energy.

Definition at line 162 of file MCGIDI_distributions.hpp.

Referenced by EnergyAngularMC().

◆ sample()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::Distributions::EnergyAngularMC::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' from the probability P(E'|E) and then samples mu from the probability P(mu|E,E'). 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 747 of file MCGIDI_headerSource.hpp.

747 {
748
749 double energyOut_1, energyOut_2;
750
751 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
752 a_input.m_energyOut1 = m_energy->sample2dOf3d( a_X, a_rng( ), a_rng, &energyOut_1, &energyOut_2 );
753 a_input.m_mu = m_angularGivenEnergy->sample( a_X, energyOut_1, energyOut_2, a_rng( ), a_rng );
754 a_input.m_phi = 2. * M_PI * a_rng( );
755 a_input.m_frame = productFrame( );
756}
#define M_PI
Definition SbMath.h:33

◆ serialize()

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

342 {
343
344 Distribution::serialize( a_buffer, a_mode );
345 m_energy = serializeProbability2d_d1( a_buffer, a_mode, m_energy );
346 m_angularGivenEnergy = serializeProbability3d( a_buffer, a_mode, m_angularGivenEnergy );
347}
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: