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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::IncoherentPhotoAtomicScattering:

Public Member Functions

LUPI_HOST_DEVICE IncoherentPhotoAtomicScattering ()
LUPI_HOST IncoherentPhotoAtomicScattering (GIDI::Distributions::IncoherentPhotoAtomicScattering const &a_incoherentPhotoAtomicScattering, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~IncoherentPhotoAtomicScattering ()
LUPI_HOST_DEVICE double energyRatio (double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double evaluateKleinNishina (double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double evaluateScatteringFactor (double a_X) 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 (Reaction const *a_reaction, LUPI_maybeUnused double a_temperature, double a_energy_in, double a_mu_lab, LUPI_maybeUnused 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 photon via incoherent photo-atomic elastic scattering.

Definition at line 278 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ IncoherentPhotoAtomicScattering() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::IncoherentPhotoAtomicScattering::IncoherentPhotoAtomicScattering ( )

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

Definition at line 801 of file MCGIDI_distributions.cc.

801 {
802
803}

◆ IncoherentPhotoAtomicScattering() [2/2]

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

811 :
812 Distribution( Type::incoherentPhotoAtomicScattering, a_incoherentPhotoAtomicScattering, a_setupInfo ) {
813
814 GUPI::Ancestry const *link = a_incoherentPhotoAtomicScattering.findInAncestry( a_incoherentPhotoAtomicScattering.href( ) );
815 GIDI::DoubleDifferentialCrossSection::IncoherentPhotoAtomicScattering const &incoherentPhotoAtomicScattering =
816 *static_cast<GIDI::DoubleDifferentialCrossSection::IncoherentPhotoAtomicScattering const *>( link );
817
818 std::string domainUnit;
819 GIDI::Functions::XYs1d const *xys1d0, *xys1d1;
820 std::size_t dataSize = 0, offset = 0;
821
822 GIDI::Functions::Function1dForm const *scatteringFactor = incoherentPhotoAtomicScattering.scatteringFactor( );
823 if( scatteringFactor->type( ) == GIDI::FormType::XYs1d ) {
824 xys1d0 = static_cast<GIDI::Functions::XYs1d const *>( scatteringFactor );
825 xys1d1 = xys1d0;
826
827 domainUnit = xys1d0->axes( )[0]->unit( );
828
829 dataSize = xys1d1->size( );
830 offset = 1; }
831 else if( scatteringFactor->type( ) == GIDI::FormType::regions1d ) {
832 GIDI::Functions::Regions1d const *regions1d = static_cast<GIDI::Functions::Regions1d const *>( scatteringFactor );
833 if( regions1d->size( ) != 2 ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor size." );
834
835 domainUnit = regions1d->axes( )[0]->unit( );
836
837 GIDI::Functions::Function1dForm const *region0 = (*regions1d)[0];
838 if( region0->type( ) != GIDI::FormType::XYs1d ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor for region 0." );
839 xys1d0 = static_cast<GIDI::Functions::XYs1d const *>( region0 );
840 if( xys1d0->size( ) != 2 ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported size of region 1 of form factor." );
841
842 GIDI::Functions::Function1dForm const *region1 = (*regions1d)[1];
843 if( region1->type( ) != GIDI::FormType::XYs1d ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor for region 1." );
844 xys1d1 = static_cast<GIDI::Functions::XYs1d const *>( region1 );
845
846 dataSize = xys1d1->size( ) + 1; }
847 else {
848 throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor. Must be XYs1d or regions1d." );
849 }
850
851 double domainFactor = 1.0;
852 if( domainUnit == "1/Ang" ) {
853 domainFactor = 0.012398419739640716; } // Converts 'h * c /Ang' to MeV.
854 else if( domainUnit == "1/cm" ) {
855 domainFactor = 0.012398419739640716 * 1e-8; } // Converts 'h * c /cm' to MeV.
856 else {
857 throw std::runtime_error( "MCGIDI::IncoherentPhotoAtomicScattering::IncoherentPhotoAtomicScattering: unsupported domain unit" );
858 }
859
860 m_energies.resize( dataSize );
861 m_scatteringFactor.resize( dataSize );
862 m_a.resize( dataSize );
863
864 std::pair<double, double> xy = (*xys1d0)[0];
865 m_energies[0] = domainFactor * xy.first;
866 m_scatteringFactor[0] = xy.second;
867 m_a[0] = 1.0;
868
869 xy = (*xys1d1)[offset];
870 double energy1 = domainFactor * xy.first;
871 double y1 = xy.second;
872
873 m_energies[1] = energy1;
874 m_scatteringFactor[1] = y1;
875
876 for( std::size_t i1 = 1 + offset; i1 < xys1d1->size( ); ++i1 ) {
877 xy = (*xys1d1)[i1];
878 double energy2 = domainFactor * xy.first;
879 double y2 = xy.second;
880
881 m_energies[i1+1-offset] = energy2;
882 m_scatteringFactor[i1+1-offset] = y2;
883
884 double _a = log( y2 / y1 ) / log( energy2 / energy1 );
885 m_a[i1-offset] = _a;
886
887 energy1 = energy2;
888 y1 = y2;
889 }
890 m_a[m_a.size()-1] = 0.0;
891}
G4ThreadLocal T * G4GeomSplitter< T >::offset
FormType type() const
Definition GIDI.hpp:667
Axes const & axes() const
Definition GIDI.hpp:1012
std::size_t size() const
Definition GIDI.hpp:1100

◆ ~IncoherentPhotoAtomicScattering()

LUPI_HOST_DEVICE MCGIDI::Distributions::IncoherentPhotoAtomicScattering::~IncoherentPhotoAtomicScattering ( )

Definition at line 896 of file MCGIDI_distributions.cc.

896 {
897
898}

Member Function Documentation

◆ angleBiasing() [1/2]

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

◆ angleBiasing() [2/2]

template<typename RNG>
LUPI_HOST_DEVICE double MCGIDI::Distributions::IncoherentPhotoAtomicScattering::angleBiasing ( Reaction const * a_reaction,
LUPI_maybeUnused double a_temperature,
double a_energy_in,
double a_mu_lab,
LUPI_maybeUnused 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 1212 of file MCGIDI_headerSource.hpp.

1213 {
1214
1215 URR_protareInfos URR_protareInfos1;
1216 double sigma = a_reaction->protareSingle( )->reactionCrossSection( a_reaction->reactionIndex( ), URR_protareInfos1, 0.0, a_energy_in );
1217
1219
1220 double one_minus_mu = 1.0 - a_mu_lab;
1221 double k_in = a_energy_in / PoPI_electronMass_MeV_c2;
1222 a_energy_out = a_energy_in / ( 1.0 + k_in * one_minus_mu );
1223 double k_out = a_energy_out / PoPI_electronMass_MeV_c2;
1224
1225 double k_ratio = k_out / k_in;
1226 double probability = evaluateScatteringFactor( a_energy_in * sqrt( 0.5 * one_minus_mu ) );
1227 probability *= k_ratio * k_ratio * ( 1.0 + a_mu_lab * a_mu_lab + k_in * k_out * one_minus_mu * one_minus_mu ) * norm;
1228
1229 return( probability );
1230}
#define MCGIDI_classicalElectronRadius
Definition MCGIDI.hpp:69
#define PoPI_electronMass_MeV_c2
Definition PoPI.hpp:31
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE double evaluateScatteringFactor(double a_X) const

◆ energyRatio()

LUPI_HOST_DEVICE double MCGIDI::Distributions::IncoherentPhotoAtomicScattering::energyRatio ( double a_energyIn,
double a_mu ) const

FIX ME.

Parameters
a_energyIn[in]
a_mu[in]

Definition at line 907 of file MCGIDI_distributions.cc.

907 {
908
909 double relativeEnergy = a_energyIn / PoPI_electronMass_MeV_c2;
910
911 return( 1.0 / ( 1.0 + relativeEnergy * ( 1.0 - a_mu ) ) );
912}

Referenced by evaluateKleinNishina().

◆ evaluateKleinNishina()

LUPI_HOST_DEVICE double MCGIDI::Distributions::IncoherentPhotoAtomicScattering::evaluateKleinNishina ( double a_energyIn,
double a_mu ) const

This method evaluates the Klein-Nishina. FIX ME. This should be a function as it does not use member data.

Parameters
a_energyIn[in]
a_mu[in]

Definition at line 921 of file MCGIDI_distributions.cc.

921 {
922
923 double relativeEnergy = a_energyIn / PoPI_electronMass_MeV_c2;
924 double _energyRatio = energyRatio( a_energyIn, a_mu );
925 double one_minus_mu = 1.0 - a_mu;
926
927 double norm = ( 1.0 + 2.0 * relativeEnergy );
928 norm = 2.0 * relativeEnergy * ( 2.0 + relativeEnergy * ( 1.0 + relativeEnergy ) * ( 8.0 + relativeEnergy ) ) / ( norm * norm );
929 norm += ( ( relativeEnergy - 2.0 ) * relativeEnergy - 2.0 ) * log( 1.0 + 2.0 * relativeEnergy );
930 norm /= relativeEnergy * relativeEnergy * relativeEnergy;
931
932 return( _energyRatio * _energyRatio * ( _energyRatio + a_mu * a_mu + relativeEnergy * one_minus_mu * one_minus_mu ) / norm );
933}
LUPI_HOST_DEVICE double energyRatio(double a_energyIn, double a_mu) const

◆ evaluateScatteringFactor()

LUPI_HOST_DEVICE double MCGIDI::Distributions::IncoherentPhotoAtomicScattering::evaluateScatteringFactor ( double a_energyIn) const

This method evaluates the Klein-Nishina.

Parameters
a_energyIn[in] The energy of the projectile.

Definition at line 941 of file MCGIDI_distributions.cc.

941 {
942
943 int intLowerIndex = binarySearchVector( a_energyIn, m_energies );
944
945 if( intLowerIndex < 1 ) {
946 if( intLowerIndex == -1 ) return( m_scatteringFactor.back( ) );
947 return( m_scatteringFactor[1] * a_energyIn / m_energies[1] );
948 }
949
950 std::size_t lowerIndex = static_cast<std::size_t>( intLowerIndex );
951 return( m_scatteringFactor[lowerIndex] * pow( a_energyIn / m_energies[lowerIndex], m_a[lowerIndex] ) );
952}
LUPI_HOST_DEVICE int binarySearchVector(double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
Definition MCGIDI.hpp:318

Referenced by angleBiasing(), and sample().

◆ sample()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::Distributions::IncoherentPhotoAtomicScattering::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 1092 of file MCGIDI_headerSource.hpp.

1092 {
1093
1094 double k1 = a_X / PoPI_electronMass_MeV_c2;
1095 double energyOut, mu, scatteringFactor;
1096
1097 if( a_X >= m_energies.back( ) ) {
1098 MCGIDI_sampleKleinNishina( k1, a_rng, &energyOut, &mu ); }
1099 else {
1100 double scatteringFactorMax = evaluateScatteringFactor( a_X );
1101 do {
1102 MCGIDI_sampleKleinNishina( k1, a_rng, &energyOut, &mu );
1103 scatteringFactor = evaluateScatteringFactor( a_X * sqrt( 0.5 * ( 1.0 - mu ) ) );
1104 } while( scatteringFactor < a_rng( ) * scatteringFactorMax );
1105 }
1106
1107 a_input.setSampledType( Sampling::SampledType::uncorrelatedBody );
1108 a_input.m_energyOut1 = energyOut * PoPI_electronMass_MeV_c2;
1109 a_input.m_mu = mu;
1110 a_input.m_phi = 2.0 * M_PI * a_rng( );
1111 a_input.m_frame = productFrame( );
1112}
LUPI_HOST_DEVICE void MCGIDI_sampleKleinNishina(double a_energyIn, RNG &&a_rng, double *a_energyOut, double *a_mu)
LUPI_HOST_DEVICE GIDI::Frame productFrame() const

◆ serialize()

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

962 {
963
964 Distribution::serialize( a_buffer, a_mode );
965
966 DATA_MEMBER_VECTOR_DOUBLE( m_energies, a_buffer, a_mode );
967 DATA_MEMBER_VECTOR_DOUBLE( m_scatteringFactor, a_buffer, a_mode );
968 DATA_MEMBER_VECTOR_DOUBLE( m_a, a_buffer, a_mode );
969}
#define DATA_MEMBER_VECTOR_DOUBLE(member, buf, mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

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