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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::Distribution:

Public Member Functions

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 is the base class for all distribution forms.

Definition at line 28 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ Distribution() [1/3]

◆ Distribution() [2/3]

LUPI_HOST MCGIDI::Distributions::Distribution::Distribution ( Type a_type,
GIDI::Distributions::Distribution const & a_distribution,
SetupInfo & a_setupInfo )
Parameters
a_type[in] The Type of the distribution.
a_distribution[in] The GIDI::Distributions::Distribution 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 43 of file MCGIDI_distributions.cc.

43 :
44 m_type( a_type ),
45 m_productFrame( a_distribution.productFrame( ) ),
46 m_projectileMass( a_setupInfo.m_protare.projectileMass( ) ),
47 m_targetMass( a_setupInfo.m_protare.targetMass( ) ),
48 m_productMass( a_setupInfo.m_product1Mass ) { // Includes nuclear excitation energy.
49
50}

◆ Distribution() [3/3]

LUPI_HOST MCGIDI::Distributions::Distribution::Distribution ( Type a_type,
GIDI::Frame a_productFrame,
SetupInfo & a_setupInfo )
Parameters
a_type[in] The Type of the distribution.
a_productFrame[in] The frame of the product's data for distribution.
a_setupInfo[in] Used internally when constructing a Protare to pass information to other constructors.

Definition at line 58 of file MCGIDI_distributions.cc.

58 :
59 m_type( a_type ),
60 m_productFrame( a_productFrame ),
61 m_projectileMass( a_setupInfo.m_protare.projectileMass( ) ),
62 m_targetMass( a_setupInfo.m_protare.targetMass( ) ),
63 m_productMass( a_setupInfo.m_product1Mass ) { // Includes nuclear excitation energy.
64
65}

◆ ~Distribution()

LUPI_HOST_DEVICE MCGIDI::Distributions::Distribution::~Distribution ( )

Definition at line 70 of file MCGIDI_distributions.cc.

70 {
71
72}

Member Function Documentation

◆ angleBiasing() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::Distribution::angleBiasing ( Reaction const * a_reaction,
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 202 of file MCGIDI_headerSource.hpp.

203 {
204
205 double probability = 0.0;
206 a_energy_out = 0.0;
207
208 switch( type( ) ) {
210 break;
212 probability = static_cast<Distributions::Unspecified const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
213 a_rng, a_energy_out );
214 break;
216 probability = static_cast<Distributions::AngularTwoBody const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
217 a_rng, a_energy_out );
218 break;
220 probability = static_cast<Distributions::KalbachMann const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
221 a_rng, a_energy_out );
222 break;
224 probability = static_cast<Distributions::Uncorrelated const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
225 a_rng, a_energy_out );
226 break;
228 probability = static_cast<Distributions::Branching3d const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
229 a_rng, a_energy_out );
230 break;
232 probability = static_cast<Distributions::EnergyAngularMC const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
233 a_rng, a_energy_out );
234 break;
236 probability = static_cast<Distributions::AngularEnergyMC const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
237 a_rng, a_energy_out );
238 break;
240 probability = static_cast<Distributions::CoherentPhotoAtomicScattering const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
241 a_rng, a_energy_out );
242 break;
244 probability = static_cast<Distributions::IncoherentPhotoAtomicScattering const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
245 a_rng, a_energy_out );
246 break;
248 probability = static_cast<Distributions::IncoherentBoundToFreePhotoAtomicScattering const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
249 a_rng, a_energy_out );
250 break;
252 probability = static_cast<Distributions::IncoherentPhotoAtomicScatteringElectron const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
253 a_rng, a_energy_out );
254 break;
256 probability = static_cast<Distributions::PairProductionGamma const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
257 a_rng, a_energy_out );
258 break;
260 probability = static_cast<Distributions::CoherentElasticTNSL const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
261 a_rng, a_energy_out );
262 break;
264 probability = static_cast<Distributions::IncoherentElasticTNSL const *>( this )->angleBiasing( a_reaction, a_temperature, a_energy_in, a_mu_lab,
265 a_rng, a_energy_out );
266 break;
267 }
268
269 return( probability );
270}
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

◆ angleBiasing() [2/2]

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

Referenced by angleBiasing().

◆ productFrame()

◆ productMass()

◆ projectileMass()

◆ sample() [1/2]

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::Distributions::Distribution::sample ( double a_X,
MCGIDI::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 135 of file MCGIDI_headerSource.hpp.

135 {
136
137 switch( type( ) ) {
139 break;
141 static_cast<Distributions::Unspecified const *>( this )->sample( a_X, a_input, a_rng );
142 break;
144 static_cast<Distributions::AngularTwoBody const *>( this )->sample( a_X, a_input, a_rng );
145 break;
147 static_cast<Distributions::KalbachMann const *>( this )->sample( a_X, a_input, a_rng );
148 break;
150 static_cast<Distributions::Uncorrelated const *>( this )->sample( a_X, a_input, a_rng );
151 break;
153 static_cast<Distributions::Branching3d const *>( this )->sample( a_X, a_input, a_rng );
154 break;
156 static_cast<Distributions::EnergyAngularMC const *>( this )->sample( a_X, a_input, a_rng );
157 break;
159 static_cast<Distributions::AngularEnergyMC const *>( this )->sample( a_X, a_input, a_rng );
160 break;
162 static_cast<Distributions::CoherentPhotoAtomicScattering const *>( this )->sample( a_X, a_input, a_rng );
163 break;
165 static_cast<Distributions::IncoherentPhotoAtomicScattering const *>( this )->sample( a_X, a_input, a_rng );
166 break;
168 static_cast<Distributions::IncoherentBoundToFreePhotoAtomicScattering const *>( this )->sample( a_X, a_input, a_rng );
169 break;
171 static_cast<Distributions::IncoherentPhotoAtomicScatteringElectron const *>( this )->sample( a_X, a_input, a_rng );
172 break;
174 static_cast<Distributions::PairProductionGamma const *>( this )->sample( a_X, a_input, a_rng );
175 break;
177 static_cast<Distributions::CoherentElasticTNSL const *>( this )->sample( a_X, a_input, a_rng );
178 break;
180 static_cast<Distributions::IncoherentElasticTNSL const *>( this )->sample( a_X, a_input, a_rng );
181 break;
182 }
183}
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL

◆ sample() [2/2]

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

Referenced by sample().

◆ serialize()

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

95 {
96
97 int distributionType = distributionTypeToInt( m_type );
98 DATA_MEMBER_INT( distributionType, a_buffer, a_mode );
99 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) m_type = intToDistributionType( distributionType );
100
101 int frame = 0;
102 if( m_productFrame == GIDI::Frame::centerOfMass ) frame = 1;
103 DATA_MEMBER_INT( frame, a_buffer, a_mode );
104 m_productFrame = GIDI::Frame::lab;
105 if( frame == 1 ) m_productFrame = GIDI::Frame::centerOfMass;
106
107 DATA_MEMBER_DOUBLE( m_projectileMass, a_buffer, a_mode );
108 DATA_MEMBER_DOUBLE( m_targetMass, a_buffer, a_mode );
109 DATA_MEMBER_DOUBLE( m_productMass, a_buffer, a_mode );
110}
#define DATA_MEMBER_DOUBLE(member, buf, mode)
#define DATA_MEMBER_INT( member, buf, mode)
@ centerOfMass
Definition GIDI.hpp:146
LUPI_HOST_DEVICE int distributionTypeToInt(Distributions::Type a_type)
LUPI_HOST_DEVICE Distributions::Type intToDistributionType(int a_type)

Referenced by MCGIDI::Distributions::AngularEnergyMC::serialize(), MCGIDI::Distributions::AngularTwoBody::serialize(), MCGIDI::Distributions::Branching3d::serialize(), MCGIDI::Distributions::CoherentElasticTNSL::serialize(), MCGIDI::Distributions::CoherentPhotoAtomicScattering::serialize(), MCGIDI::Distributions::EnergyAngularMC::serialize(), MCGIDI::Distributions::IncoherentBoundToFreePhotoAtomicScattering::serialize(), MCGIDI::Distributions::IncoherentElasticTNSL::serialize(), MCGIDI::Distributions::IncoherentPhotoAtomicScattering::serialize(), MCGIDI::Distributions::IncoherentPhotoAtomicScatteringElectron::serialize(), MCGIDI::Distributions::KalbachMann::serialize(), MCGIDI::Distributions::PairProductionGamma::serialize(), MCGIDI::Distributions::Uncorrelated::serialize(), and MCGIDI::Distributions::Unspecified::serialize().

◆ setModelDBRC_data()

LUPI_HOST void MCGIDI::Distributions::Distribution::setModelDBRC_data ( Sampling::Upscatter::ModelDBRC_data * a_modelDBRC_data)

This method calls the **setModelDBRC_data2* method if the distribution is AngularTwoBody, otherwise it * executes a throw.

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

Definition at line 80 of file MCGIDI_distributions.cc.

80 {
81
82 if( type( ) != Type::angularTwoBody ) throw std::runtime_error( "Setting ModelDBRC_data for non-two body distribution is not allowed." );
83
84 static_cast<AngularTwoBody *>( this )->setModelDBRC_data2( a_modelDBRC_data );
85}

◆ targetMass()

◆ type()

LUPI_HOST_DEVICE Type MCGIDI::Distributions::Distribution::type ( ) const
inline

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