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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::KalbachMann:

Public Member Functions

LUPI_HOST_DEVICE KalbachMann ()
LUPI_HOST KalbachMann (GIDI::Distributions::KalbachMann const &a_KalbachMann, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~KalbachMann ()
LUPI_HOST_DEVICE double energyToMeVFactor () const
LUPI_HOST_DEVICE double eb_massFactor () const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1f () const
LUPI_HOST_DEVICE Functions::Function2dr () const
LUPI_HOST_DEVICE Functions::Function2da () 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)
LUPI_HOST_DEVICE double evaluate (double E_in_lab, double E_out, double mu)
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 whose distribution is represented by Kalbach-Mann systematics.

Definition at line 203 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ KalbachMann() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::KalbachMann::KalbachMann ( )

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

Definition at line 410 of file MCGIDI_distributions.cc.

410 :
411 m_energyToMeVFactor( 0.0 ),
412 m_eb_massFactor( 0.0 ),
413 m_f( nullptr ),
414 m_r( nullptr ),
415 m_a( nullptr ) {
416
417}

Referenced by KalbachMann().

◆ KalbachMann() [2/2]

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

424 :
425 Distribution( Type::KalbachMann, a_KalbachMann, a_setupInfo ),
426 m_energyToMeVFactor( 1 ), // FIXME.
427 m_eb_massFactor( 1 ), // FIXME.
428 m_f( Probabilities::parseProbability2d_d1( a_KalbachMann.f( ), nullptr ) ),
429 m_r( Functions::parseFunction2d( a_KalbachMann.r( ) ) ),
430 m_a( Functions::parseFunction2d( a_KalbachMann.a( ) ) ) {
431
432}
LUPI_HOST Function2d * parseFunction2d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST ProbabilityBase2d_d1 * parseProbability2d_d1(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)

◆ ~KalbachMann()

LUPI_HOST_DEVICE MCGIDI::Distributions::KalbachMann::~KalbachMann ( )

Definition at line 437 of file MCGIDI_distributions.cc.

437 {
438
439 delete m_f;
440 delete m_r;
441 delete m_a;
442}

Member Function Documentation

◆ a()

LUPI_HOST_DEVICE Functions::Function2d * MCGIDI::Distributions::KalbachMann::a ( ) const
inline

Returns the value of the m_a.

Definition at line 221 of file MCGIDI_distributions.hpp.

Referenced by KalbachMann().

◆ angleBiasing() [1/2]

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

917 {
918
919 a_energy_out = 0.0;
920
921 double initialMass = projectileMass( ) + targetMass( );
922 double energy_out_com = m_f->sample( a_energy_in, a_rng( ), a_rng );
923 double productBeta = MCGIDI_particleBeta( productMass( ), energy_out_com );
924 double boostBeta = sqrt( a_energy_in * ( a_energy_in + 2. * projectileMass( ) ) ) / ( a_energy_in + initialMass ); // beta = v/c.
925
926 double muPlus = 0.0, JacobianPlus = 0.0, muMinus = 0.0, JacobianMinus = 0.0;
927
928 int numberOfMus = muCOM_From_muLab( a_mu_lab, boostBeta, productBeta, muPlus, JacobianPlus, muMinus, JacobianMinus );
929
930 if( numberOfMus == 0 ) return( 0.0 );
931
932 double rAtEnergyEnergyPrime = m_r->evaluate( a_energy_in, energy_out_com );
933 double aAtEnergyEnergyPrime = m_a->evaluate( a_energy_in, energy_out_com );
934 double aMu = aAtEnergyEnergyPrime * muPlus;
935
936 double probability = 0.5 * JacobianPlus;
937 if( productMass( ) == 0.0 ) {
938 probability *= 1.0 - rAtEnergyEnergyPrime + rAtEnergyEnergyPrime * aAtEnergyEnergyPrime * exp( aMu ) / sinh( aAtEnergyEnergyPrime ); }
939 else {
940 probability *= aAtEnergyEnergyPrime * ( cosh( aMu ) + rAtEnergyEnergyPrime * cosh( aMu ) ) / sinh( aAtEnergyEnergyPrime );
941 }
942
943 if( numberOfMus == 2 ) {
944 aMu = aAtEnergyEnergyPrime * muMinus;
945
946 double probabilityMinus = 0.5 * JacobianMinus;
947 if( productMass( ) == 0.0 ) {
948 probabilityMinus *= 1.0 - rAtEnergyEnergyPrime + rAtEnergyEnergyPrime * aAtEnergyEnergyPrime * exp( aMu ) / sinh( aAtEnergyEnergyPrime ); }
949 else {
950 probabilityMinus *= aAtEnergyEnergyPrime * ( cosh( aMu ) + rAtEnergyEnergyPrime * cosh( aMu ) ) / sinh( aAtEnergyEnergyPrime );
951 }
952 probability += probabilityMinus;
953
954 if( probabilityMinus > a_rng( ) * probability ) muPlus = muMinus;
955 }
956
957 double productBeta2 = productBeta * productBeta;
958 double productBetaLab2 = productBeta2 + boostBeta * boostBeta * ( 1.0 - productBeta2 * ( 1.0 - muPlus * muPlus ) ) + 2.0 * muPlus * productBeta * boostBeta;
959 productBetaLab2 /= 1.0 - muPlus * productBeta * boostBeta;
960 a_energy_out = MCGIDI::particleKineticEnergyFromBeta2( productMass( ), productBetaLab2 );
961
962 return( probability );
963}
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71
LUPI_HOST_DEVICE double targetMass() 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::KalbachMann::angleBiasing ( Reaction const * a_reaction,
double a_temperature,
double a_energy_in,
double a_mu_lab,
RNG && a_rng,
double & a_energy_out ) const

◆ eb_massFactor()

LUPI_HOST_DEVICE double MCGIDI::Distributions::KalbachMann::eb_massFactor ( ) const
inline

Returns the value of the m_eb_massFactor.

Definition at line 218 of file MCGIDI_distributions.hpp.

◆ energyToMeVFactor()

LUPI_HOST_DEVICE double MCGIDI::Distributions::KalbachMann::energyToMeVFactor ( ) const
inline

Returns the value of the m_energyToMeVFactor.

Definition at line 217 of file MCGIDI_distributions.hpp.

◆ evaluate()

LUPI_HOST_DEVICE double MCGIDI::Distributions::KalbachMann::evaluate ( double a_energy,
double a_energyOut,
double a_mu )

This method evaluates the Kalbach-Mann formalism at the projectile energy a_energy, and outgoing product energy a_energyOut and a_mu.

Parameters
a_energy[in] The energy of the projectile in the lab frame.
a_energyOut[in] The energy of the product in the center-of-mass frame.
a_mu[in] The mu of the product in the center-of-mass frame.

Definition at line 452 of file MCGIDI_distributions.cc.

452 {
453
454// double f_0 = m_f->evaluate( a_energy, a_energyOut );
455 double rValue = m_r->evaluate( a_energy, a_energyOut );
456 double aValue = m_a->evaluate( a_energy, a_energyOut );
457// double pdf_val = aValue * f_0 / 2.0 / sinh( aValue ) * (cosh(aValue * a_mu) + rValue * sinh( aValue * a_mu ) ); // double-differential PDF for a_energyOut and a_mu (Eq. 6.4 in ENDF-102, 2012)
458 double pdf_val = aValue / ( 2.0 * sinh( aValue ) ) * ( cosh( aValue * a_mu ) + rValue * sinh( aValue * a_mu ) ); // double-differential PDF for a_energyOut and mu (Eq. 6.4 in ENDF-102, 2012)
459 return pdf_val;
460}

◆ f()

LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * MCGIDI::Distributions::KalbachMann::f ( ) const
inline

Returns the value of the m_f.

Definition at line 219 of file MCGIDI_distributions.hpp.

Referenced by KalbachMann().

◆ r()

LUPI_HOST_DEVICE Functions::Function2d * MCGIDI::Distributions::KalbachMann::r ( ) const
inline

Returns the value of the m_r.

Definition at line 220 of file MCGIDI_distributions.hpp.

Referenced by KalbachMann().

◆ sample()

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

This method samples the outgoing product data using the Kalbach-Mann formalism.

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 877 of file MCGIDI_headerSource.hpp.

877 {
878
879 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
880 a_input.m_energyOut1 = m_f->sample( a_X, a_rng( ), a_rng );
881 double rValue = m_r->evaluate( a_X, a_input.m_energyOut1 );
882 double aValue = m_a->evaluate( a_X, a_input.m_energyOut1 );
883
884 // In the following: Cosh[ a mu ] + r Sinh[ a mu ] = ( 1 - r ) Cosh[ a mu ] + r ( Cosh[ a mu ] + Sinh[ a mu ] ).
885 if( a_rng( ) >= rValue ) { // Sample the '( 1 - r ) Cosh[ a mu ]' term.
886 double T = ( 2. * a_rng( ) - 1. ) * sinh( aValue );
887
888 a_input.m_mu = log( T + sqrt( T * T + 1. ) ) / aValue; }
889 else { // Sample the 'r ( Cosh[ a mu ] + Sinh[ a mu ] )' term.
890 double rng1 = a_rng( ), exp_a = exp( aValue );
891
892 a_input.m_mu = log( rng1 * exp_a + ( 1. - rng1 ) / exp_a ) / aValue;
893 }
894 if( a_input.m_mu < -1 ) a_input.m_mu = -1;
895 if( a_input.m_mu > 1 ) a_input.m_mu = 1;
896
897 a_input.m_phi = 2. * M_PI * a_rng( );
898 a_input.m_frame = productFrame( );
899}
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE GIDI::Frame productFrame() const

◆ serialize()

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

470 {
471
472 Distribution::serialize( a_buffer, a_mode );
473 DATA_MEMBER_DOUBLE( m_energyToMeVFactor, a_buffer, a_mode );
474 DATA_MEMBER_DOUBLE( m_eb_massFactor, a_buffer, a_mode );
475
476 m_f = serializeProbability2d_d1( a_buffer, a_mode, m_f );
477 m_r = serializeFunction2d( a_buffer, a_mode, m_r );
478 m_a = serializeFunction2d( a_buffer, a_mode, m_a );
479}
#define DATA_MEMBER_DOUBLE(member, buf, mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Functions::Function2d * serializeFunction2d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function2d *a_function2d)
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: