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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::AngularTwoBody:

Public Member Functions

LUPI_HOST_DEVICE AngularTwoBody ()
LUPI_HOST AngularTwoBody (GIDI::Distributions::AngularTwoBody const &a_angularTwoBody, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~AngularTwoBody ()
LUPI_HOST_DEVICE double residualMass () const
LUPI_HOST_DEVICE double Q () const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1angular () 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 bool Upscatter () const
LUPI_HOST void setModelDBRC_data2 (Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
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 a two-body interaction.

Definition at line 66 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ AngularTwoBody() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::AngularTwoBody::AngularTwoBody ( )

Base contructor.

Definition at line 120 of file MCGIDI_distributions.cc.

120 :
121 m_residualMass( 0.0 ),
122 m_Q( 0.0 ),
123 m_twoBodyThreshold( 0.0 ),
124 m_Upscatter( false ),
125 m_angular( nullptr ),
126 m_modelDBRC_data( nullptr ) {
127
128}

◆ AngularTwoBody() [2/2]

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

135 :
136 Distribution( Type::angularTwoBody, a_angularTwoBody, a_setupInfo ),
137 m_residualMass( a_setupInfo.m_product2Mass ), // Includes nuclear excitation energy.
138 m_Q( a_setupInfo.m_Q ),
139 m_twoBodyThreshold( a_setupInfo.m_reaction->twoBodyThreshold( ) ),
140 m_Upscatter( false ),
141 m_angular( Probabilities::parseProbability2d_d1( a_angularTwoBody.angular( ), &a_setupInfo ) ),
142 m_modelDBRC_data( nullptr ) {
143
144 if( a_setupInfo.m_protare.projectileIntid( ) == PoPI::Intids::neutron ) {
145 m_Upscatter = a_setupInfo.m_reaction->ENDF_MT( ) == 2;
146 }
147}
LUPI_HOST ProbabilityBase2d_d1 * parseProbability2d_d1(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)
static int constexpr neutron
Definition PoPI.hpp:177

◆ ~AngularTwoBody()

LUPI_HOST_DEVICE MCGIDI::Distributions::AngularTwoBody::~AngularTwoBody ( )

Definition at line 152 of file MCGIDI_distributions.cc.

152 {
153
154 delete m_angular;
155 delete m_modelDBRC_data;
156}

Member Function Documentation

◆ angleBiasing() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::AngularTwoBody::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[out] The energy of the emitted outgoing particle.
Returns
The probability of emitting outgoing particle into lab angle a_mu_lab.

Definition at line 410 of file MCGIDI_headerSource.hpp.

411 {
412
413 a_energy_out = 0.0;
414
415 double initialMass = projectileMass( ) + targetMass( ), finalMass = productMass( ) + m_residualMass;
416 double _x = targetMass( ) * ( a_energy_in - m_twoBodyThreshold ) / ( finalMass * finalMass );
417 double Kp; // Total kinetic energy of products in the center-of-mass.
418
419 if( _x < 2e-5 ) {
420 Kp = finalMass * _x * ( 1 - 0.5 * _x * ( 1 - _x ) ); }
421 else { // This is the relativistic formula derived from E^2 - (pc)^2 which is frame independent (i.e., an invariant).
422 Kp = sqrt( finalMass * finalMass + 2.0 * targetMass( ) * ( a_energy_in - m_twoBodyThreshold ) ) - finalMass;
423 }
424 if( Kp < 0 ) Kp = 0.; // FIXME There needs to be a better test here.
425
426 double energy_product_com = 0.5 * Kp * ( Kp + 2.0 * m_residualMass ) / ( Kp + productMass( ) + m_residualMass );
427
428 if( productMass( ) == 0.0 ) {
429 double boostBeta = sqrt( a_energy_in * ( a_energy_in + 2. * projectileMass( ) ) ) / ( a_energy_in + initialMass ); // Good, even for projectileMass = 0.
430 double one_mu_beta = 1.0 - a_mu_lab * boostBeta;
431 double mu_com = ( a_mu_lab - boostBeta ) / one_mu_beta;
432 double Jacobian = ( 1.0 - boostBeta * boostBeta ) / ( one_mu_beta * one_mu_beta );
433
434 a_energy_out = sqrt( 1.0 - boostBeta * boostBeta ) * energy_product_com * ( 1.0 + mu_com * boostBeta );
435
436 return( Jacobian * m_angular->evaluate( a_energy_in, mu_com ) );
437 }
438
439 double productBeta = MCGIDI_particleBeta( productMass( ), energy_product_com );
440 double boostBeta = sqrt( a_energy_in * ( a_energy_in + 2. * projectileMass( ) ) ) / ( a_energy_in + initialMass ); // beta = v/c.
441 double muPlus = 0.0, JacobianPlus = 0.0, muMinus = 0.0, JacobianMinus = 0.0;
442
443 int numberOfMus = muCOM_From_muLab( a_mu_lab, boostBeta, productBeta, muPlus, JacobianPlus, muMinus, JacobianMinus );
444
445 if( numberOfMus == 0 ) return( 0.0 );
446
447 double probability = JacobianPlus * m_angular->evaluate( a_energy_in, muPlus );
448
449 if( numberOfMus == 2 ) {
450 double probabilityMinus = JacobianMinus * m_angular->evaluate( a_energy_in, muMinus );
451 probability += probabilityMinus;
452 if( probabilityMinus > a_rng( ) * probability ) {
453 muPlus = muMinus;
454 }
455 }
456
457 double productBeta2 = productBeta * productBeta;
458 double productBetaLab2 = productBeta2 + boostBeta * boostBeta * ( 1.0 - productBeta2 * ( 1.0 - muPlus * muPlus ) ) + 2.0 * muPlus * productBeta * boostBeta;
459 productBetaLab2 /= 1.0 - muPlus * productBeta * boostBeta;
460 a_energy_out = MCGIDI::particleKineticEnergyFromBeta2( productMass( ), productBetaLab2 );
461
462 return( probability );
463}
#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::AngularTwoBody::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::AngularTwoBody::angular ( ) const
inline

Returns the value of the m_angular.

Definition at line 86 of file MCGIDI_distributions.hpp.

Referenced by AngularTwoBody().

◆ Q()

LUPI_HOST_DEVICE double MCGIDI::Distributions::AngularTwoBody::Q ( ) const
inline

Returns the value of the m_Q.

Definition at line 85 of file MCGIDI_distributions.hpp.

◆ residualMass()

LUPI_HOST_DEVICE double MCGIDI::Distributions::AngularTwoBody::residualMass ( ) const
inline

Returns the value of the m_residualMass.

Definition at line 84 of file MCGIDI_distributions.hpp.

◆ sample()

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

This method samples the outgoing product data for the two outgoing particles in a two-body outgoing channel. First, is samples mu, the cosine of the product's outgoing angle, since this is for two-body interactions, mu is in the center-of-mass frame. It then calls kinetics_COMKineticEnergy2LabEnergyAndMomentum.

Parameters
a_X[in] The energy of the projectile in the lab frame.
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 366 of file MCGIDI_headerSource.hpp.

366 {
367
368 double initialMass = projectileMass( ) + targetMass( ), finalMass = productMass( ) + m_residualMass;
369 double beta = sqrt( a_X * ( a_X + 2. * projectileMass( ) ) ) / ( a_X + initialMass ); // beta = v/c.
370 double _x = targetMass( ) * ( a_X - m_twoBodyThreshold ) / ( finalMass * finalMass );
371 double Kp; // Kp is the total kinetic energy for m3 and m4 in the COM frame.
372
373 a_input.setSampledType( Sampling::SampledType::firstTwoBody );
374
375 if( m_Upscatter ) {
376 if( ( a_input.m_upscatterModel == Sampling::Upscatter::Model::B ) || ( a_input.m_upscatterModel == Sampling::Upscatter::Model::BSnLimits )
377 || ( a_input.m_upscatterModel == Sampling::Upscatter::Model::DBRC ) ) {
378 if( upscatterModelB( a_X, a_input, a_rng ) ) return;
379 }
380 }
381
382 if( _x < 2e-5 ) {
383 Kp = finalMass * _x * ( 1 - 0.5 * _x * ( 1 - _x ) ); }
384 else { // This is the relativistic formula derived from E^2 - (pc)^2 is frame independent.
385 Kp = sqrt( finalMass * finalMass + 2 * targetMass( ) * ( a_X - m_twoBodyThreshold ) ) - finalMass;
386 }
387 if( Kp < 0 ) Kp = 0.; // FIXME There needs to be a better test here.
388
389 a_input.m_mu = m_angular->sample( a_X, a_rng( ), a_rng );
390 a_input.m_phi = 2. * M_PI * a_rng( );
391 kinetics_COMKineticEnergy2LabEnergyAndMomentum( beta, Kp, productMass( ), m_residualMass, a_input );
392}
LUPI_HOST_DEVICE void kinetics_COMKineticEnergy2LabEnergyAndMomentum(double a_beta, double a_kinetic_com, double a_m3cc, double a_m4cc, MCGIDI::Sampling::Input &a_input)
#define M_PI
Definition SbMath.h:33

◆ serialize()

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

166 {
167
168 Distribution::serialize( a_buffer, a_mode );
169 DATA_MEMBER_DOUBLE( m_residualMass, a_buffer, a_mode );
170 DATA_MEMBER_DOUBLE( m_Q, a_buffer, a_mode );
171 DATA_MEMBER_DOUBLE( m_twoBodyThreshold, a_buffer, a_mode );
172 DATA_MEMBER_INT( m_Upscatter, a_buffer, a_mode );
173
174 m_angular = serializeProbability2d_d1( a_buffer, a_mode, m_angular );
175 m_modelDBRC_data = serializeModelDBRC_data( a_buffer, a_mode, m_modelDBRC_data );
176}
#define DATA_MEMBER_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)
LUPI_HOST_DEVICE ModelDBRC_data * serializeModelDBRC_data(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, ModelDBRC_data *a_modelDBRC_data)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * serializeProbability2d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d1 *a_probability2d)

◆ setModelDBRC_data2()

LUPI_HOST void MCGIDI::Distributions::AngularTwoBody::setModelDBRC_data2 ( Sampling::Upscatter::ModelDBRC_data * a_modelDBRC_data)

This method sets this m_modelDBRC_data to a_modelDBRC_data. It also deletes the current m_modelDBRC_data member.

Parameters
a_modelDBRC_data[in] The instance storing data needed to treat the DRRC upscatter mode.

Definition at line 184 of file MCGIDI_distributions.cc.

184 {
185
186 delete m_modelDBRC_data;
187 m_modelDBRC_data = a_modelDBRC_data;
188}

◆ Upscatter()

LUPI_HOST_DEVICE bool MCGIDI::Distributions::AngularTwoBody::Upscatter ( ) const
inline

Returns the value of the m_Upscatter.

Definition at line 93 of file MCGIDI_distributions.hpp.


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