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

#include <MCGIDI_distributions.hpp>

Inheritance diagram for MCGIDI::Distributions::CoherentPhotoAtomicScattering:

Public Member Functions

LUPI_HOST_DEVICE CoherentPhotoAtomicScattering ()
LUPI_HOST CoherentPhotoAtomicScattering (GIDI::Distributions::CoherentPhotoAtomicScattering const &a_coherentPhotoAtomicScattering, SetupInfo &a_setupInfo)
LUPI_HOST_DEVICE ~CoherentPhotoAtomicScattering ()
LUPI_HOST_DEVICE double evaluate (double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double evaluateFormFactor (double a_energyIn, double a_mu) 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 coherent photo-atomic elastic scattering.

Definition at line 237 of file MCGIDI_distributions.hpp.

Constructor & Destructor Documentation

◆ CoherentPhotoAtomicScattering() [1/2]

LUPI_HOST_DEVICE MCGIDI::Distributions::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering ( )

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

Definition at line 489 of file MCGIDI_distributions.cc.

489 :
490 m_realAnomalousFactor( nullptr ),
491 m_imaginaryAnomalousFactor( nullptr ) {
492
493}

◆ CoherentPhotoAtomicScattering() [2/2]

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

500 :
501 Distribution( Type::coherentPhotoAtomicScattering, a_coherentPhotoAtomicScattering, a_setupInfo ),
502 m_anomalousDataPresent( false ),
503 m_realAnomalousFactor( nullptr ),
504 m_imaginaryAnomalousFactor( nullptr ) {
505
506 GUPI::Ancestry const *link = a_coherentPhotoAtomicScattering.findInAncestry( a_coherentPhotoAtomicScattering.href( ) );
507 GIDI::DoubleDifferentialCrossSection::CoherentPhotoAtomicScattering const &coherentPhotoAtomicScattering =
508 *static_cast<GIDI::DoubleDifferentialCrossSection::CoherentPhotoAtomicScattering const *>( link );
509
510 std::string domainUnit;
511 GIDI::Functions::XYs1d const *xys1d0, *xys1d1;
512 std::size_t dataSize = 0, offset = 0;
513
514 GIDI::Functions::Function1dForm const *formFactor = coherentPhotoAtomicScattering.formFactor( );
515 if( formFactor->type( ) == GIDI::FormType::XYs1d ) {
516 xys1d0 = static_cast<GIDI::Functions::XYs1d const *>( formFactor );
517 xys1d1 = xys1d0;
518
519 domainUnit = xys1d0->axes( )[0]->unit( );
520
521 dataSize = xys1d1->size( );
522 offset = 1; }
523 else if( formFactor->type( ) == GIDI::FormType::regions1d ) {
524 GIDI::Functions::Regions1d const *regions1d = static_cast<GIDI::Functions::Regions1d const *>( formFactor );
525 if( regions1d->size( ) != 2 ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor size." );
526
527 domainUnit = regions1d->axes( )[0]->unit( );
528
529 GIDI::Functions::Function1dForm const *region0 = (*regions1d)[0];
530 if( region0->type( ) != GIDI::FormType::XYs1d ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor for region 0." );
531 xys1d0 = static_cast<GIDI::Functions::XYs1d const *>( region0 );
532 if( xys1d0->size( ) != 2 ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported size of region 1 of form factor." );
533
534 GIDI::Functions::Function1dForm const *region1 = (*regions1d)[1];
535 if( region1->type( ) != GIDI::FormType::XYs1d ) throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor for region 1." );
536 xys1d1 = static_cast<GIDI::Functions::XYs1d const *>( region1 );
537
538 dataSize = xys1d1->size( ) + 1; }
539 else {
540 throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported form factor. Must be XYs1d or regions1d." );
541 }
542
543 double domainFactor = 1.0;
544 if( domainUnit == "1/Ang" ) {
545 domainFactor = 0.012398419739640716; } // Converts 'h * c /Ang' to MeV.
546 else if( domainUnit == "1/cm" ) {
547 domainFactor = 0.012398419739640716 * 1e-8; } // Converts 'h * c /cm' to MeV.
548 else {
549 throw std::runtime_error( "MCGIDI::CoherentPhotoAtomicScattering::CoherentPhotoAtomicScattering: unsupported domain unit" );
550 }
551
552 m_energies.resize( dataSize );
553 m_formFactor.resize( dataSize );
554 m_a.resize( dataSize );
555 m_integratedFormFactor.resize( dataSize );
556 m_integratedFormFactorSquared.resize( dataSize );
557 m_probabilityNorm1_1.resize( dataSize );
558 m_probabilityNorm1_3.resize( dataSize );
559 m_probabilityNorm1_5.resize( dataSize );
560 m_probabilityNorm2_1.resize( dataSize );
561 m_probabilityNorm2_3.resize( dataSize );
562 m_probabilityNorm2_5.resize( dataSize );
563
564 std::pair<double, double> xy = (*xys1d0)[0];
565 m_energies[0] = 0.0;
566 m_formFactor[0] = xy.second;
567 m_a[0] = 0.0;
568 m_integratedFormFactor[0] = 0.0;
569 m_integratedFormFactorSquared[0] = 0.0;
570
571 xy = (*xys1d1)[offset];
572 double energy1 = domainFactor * xy.first;
573 double y1 = xy.second;
574 m_energies[1] = energy1;
575 m_formFactor[1] = y1;
576 m_integratedFormFactor[1] = 0.5 * energy1 * energy1 * y1;
577 m_integratedFormFactorSquared[1] = 0.5 * energy1 * energy1 * y1 * y1;
578
579 double sum1 = m_integratedFormFactor[1];
580 double sum2 = m_integratedFormFactorSquared[1];
581 for( std::size_t i1 = 1 + offset; i1 < xys1d1->size( ); ++i1 ) {
582 xy = (*xys1d1)[i1];
583 double energy2 = domainFactor * xy.first;
584 double y2 = xy.second;
585
586 double logEs = log( energy2 / energy1 );
587 double _a = log( y2 / y1 ) / logEs;
588
589 m_energies[i1+1-offset] = energy2;
590 m_formFactor[i1+1-offset] = y2;
591 m_a[i1-offset] = _a;
592
593 sum1 += coherentPhotoAtomicScatteringIntegrateSub( 1, _a, logEs, energy1, y1, energy2, y2 );
594 m_integratedFormFactor[i1+1-offset] = sum1;
595
596 sum2 += coherentPhotoAtomicScatteringIntegrateSub( 1, 2.0 * _a, logEs, energy1, y1 * y1, energy2, y2 * y2 );
597 m_integratedFormFactorSquared[i1+1-offset] = sum2;
598
599 energy1 = energy2;
600 y1 = y2;
601 }
602
603 m_a[m_a.size()-1] = 0.0;
604
605 if( coherentPhotoAtomicScattering.realAnomalousFactor( ) != nullptr ) {
606 m_anomalousDataPresent = true;
607 m_realAnomalousFactor = Functions::parseFunction1d_d1( coherentPhotoAtomicScattering.realAnomalousFactor( ) );
608 m_imaginaryAnomalousFactor = Functions::parseFunction1d_d1( coherentPhotoAtomicScattering.imaginaryAnomalousFactor( ) );
609 }
610
611 m_probabilityNorm1_1[0] = 0.0;
612 m_probabilityNorm1_3[0] = 0.0;
613 m_probabilityNorm1_5[0] = 0.0;
614 m_probabilityNorm2_1[0] = 0.0;
615 m_probabilityNorm2_3[0] = 0.0;
616 m_probabilityNorm2_5[0] = 0.0;
617 energy1 = m_energies[1];
618 y1 = m_formFactor[0];
619 for( std::size_t i1 = 1; i1 < m_probabilityNorm1_1.size( ); ++i1 ) {
620 double energy2 = m_energies[i1];
621 double y2 = m_formFactor[i1];
622 double logEs = log( energy2 / energy1 );
623
624 m_probabilityNorm1_1[i1] = m_probabilityNorm1_1[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 1, m_a[i1-1], logEs, energy1, y1, energy2, y2 );
625 m_probabilityNorm1_3[i1] = m_probabilityNorm1_3[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 3, m_a[i1-1], logEs, energy1, y1, energy2, y2 );
626 m_probabilityNorm1_5[i1] = m_probabilityNorm1_5[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 5, m_a[i1-1], logEs, energy1, y1, energy2, y2 );
627
628 m_probabilityNorm2_1[i1] = m_probabilityNorm2_1[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 1, 2.0 * m_a[i1-1], logEs, energy1, y1 * y1, energy2, y2 * y2 );
629 m_probabilityNorm2_3[i1] = m_probabilityNorm2_3[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 3, 2.0 * m_a[i1-1], logEs, energy1, y1 * y1, energy2, y2 * y2 );
630 m_probabilityNorm2_5[i1] = m_probabilityNorm2_5[i1-1] + coherentPhotoAtomicScatteringIntegrateSub( 5, 2.0 * m_a[i1-1], logEs, energy1, y1 * y1, energy2, y2 * y2 );
631
632 energy1 = energy2;
633 y1 = y2;
634 }
635}
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
LUPI_HOST Function1d_d1 * parseFunction1d_d1(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)

◆ ~CoherentPhotoAtomicScattering()

LUPI_HOST_DEVICE MCGIDI::Distributions::CoherentPhotoAtomicScattering::~CoherentPhotoAtomicScattering ( )

Definition at line 640 of file MCGIDI_distributions.cc.

640 {
641
642 delete m_realAnomalousFactor;
643 delete m_imaginaryAnomalousFactor;
644}

Member Function Documentation

◆ angleBiasing() [1/2]

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

1062 {
1063
1064 a_energy_out = a_energy_in;
1065
1066 URR_protareInfos URR_protareInfos1;
1067 double sigma = a_reaction->protareSingle( )->reactionCrossSection( a_reaction->reactionIndex( ), URR_protareInfos1, 0.0, a_energy_in );
1068 double formFactor = evaluateFormFactor( a_energy_in, a_mu_lab );
1069 double imaginaryAnomalousFactor = 0.0;
1070
1071 if( m_anomalousDataPresent ) {
1072 formFactor += m_realAnomalousFactor->evaluate( a_energy_in );
1073 imaginaryAnomalousFactor = m_imaginaryAnomalousFactor->evaluate( a_energy_in );
1074 }
1075
1076 double probability = M_PI * MCGIDI_classicalElectronRadius * MCGIDI_classicalElectronRadius * ( 1.0 + a_mu_lab * a_mu_lab )
1077 * ( formFactor * formFactor + imaginaryAnomalousFactor * imaginaryAnomalousFactor ) / sigma;
1078
1079 return( probability );
1080}
#define MCGIDI_classicalElectronRadius
Definition MCGIDI.hpp:69
#define M_PI
Definition SbMath.h:33
LUPI_HOST_DEVICE double evaluateFormFactor(double a_energyIn, double a_mu) const

◆ evaluate()

LUPI_HOST_DEVICE double MCGIDI::Distributions::CoherentPhotoAtomicScattering::evaluate ( double a_energyIn,
double a_mu ) const

This method evaluates the coherent photo-atomic scattering double differentil at the projectile energy a_energy and product cosine of angle a_mu.

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

Definition at line 653 of file MCGIDI_distributions.cc.

653 {
654
655 double probability;
656 int intLowerIndexEnergy = binarySearchVector( a_energyIn, m_energies, true ); // FIXME - need to handle case where lowerIndexEnergy = 0 like in evaluateScatteringFactor.
657 std::size_t lowerIndexEnergy = static_cast<std::size_t>( intLowerIndexEnergy );
658 double _a = m_a[lowerIndexEnergy];
659 double _a_2 = _a * _a;
660 double X1 = m_energies[lowerIndexEnergy];
661 double logEs = log( a_energyIn / X1 );
662 double formFactor_1 = m_formFactor[lowerIndexEnergy];
663 double formFactor_2 = formFactor_1 * formFactor_1;
664 double formFactorEnergyIn_1 = formFactor_1 * pow( a_energyIn / X1, _a );
665 double formFactorEnergyIn_2 = formFactorEnergyIn_1 * formFactorEnergyIn_1;
666 double inverseEnergyIn_1 = 1.0 / a_energyIn;
667 double inverseEnergyIn_2 = inverseEnergyIn_1 * inverseEnergyIn_1;
668 double inverseEnergyIn_3 = inverseEnergyIn_1 * inverseEnergyIn_2;
669 double inverseEnergyIn_4 = inverseEnergyIn_2 * inverseEnergyIn_2;
670 double inverseEnergyIn_5 = inverseEnergyIn_1 * inverseEnergyIn_4;
671 double inverseEnergyIn_6 = inverseEnergyIn_2 * inverseEnergyIn_4;
672
673 double norm = 0.5 * inverseEnergyIn_2 * ( m_probabilityNorm2_1[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 1, _a_2, logEs, X1, formFactor_2, a_energyIn, formFactorEnergyIn_2 ) )
674 - inverseEnergyIn_4 * ( m_probabilityNorm2_3[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 3, _a_2, logEs, X1, formFactor_2, a_energyIn, formFactorEnergyIn_2 ) )
675 + inverseEnergyIn_6 * ( m_probabilityNorm2_5[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 5, _a_2, logEs, X1, formFactor_2, a_energyIn, formFactorEnergyIn_2 ) );
676
677 double realAnomalousFactor = 0.0;
678 double imaginaryAnomalousFactor = 0.0;
679 if( m_anomalousDataPresent ) {
680 realAnomalousFactor = m_realAnomalousFactor->evaluate( a_energyIn );
681 imaginaryAnomalousFactor = m_imaginaryAnomalousFactor->evaluate( a_energyIn );
682 norm += realAnomalousFactor * ( inverseEnergyIn_1 * ( m_probabilityNorm1_1[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 1, _a, logEs, X1, formFactor_1, a_energyIn, formFactorEnergyIn_1 ) )
683 - 2.0 * inverseEnergyIn_3 * ( m_probabilityNorm1_3[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 3, _a, logEs, X1, formFactor_1, a_energyIn, formFactorEnergyIn_1 ) )
684 + 2.0 * inverseEnergyIn_5 * ( m_probabilityNorm1_5[lowerIndexEnergy] + coherentPhotoAtomicScatteringIntegrateSub( 5, _a, logEs, X1, formFactor_1, a_energyIn, formFactorEnergyIn_1 ) ) );
685 }
686 norm *= 16.0;
687 norm += 8.0 / 3.0 * ( realAnomalousFactor * realAnomalousFactor + imaginaryAnomalousFactor * imaginaryAnomalousFactor );
688
689 double _formFactor = evaluateFormFactor( a_energyIn, a_mu );
690 probability = ( 1.0 + a_mu * a_mu ) * ( ( _formFactor + realAnomalousFactor ) * ( _formFactor + realAnomalousFactor ) + imaginaryAnomalousFactor * imaginaryAnomalousFactor ) / norm;
691
692 return( probability );
693}
LUPI_HOST_DEVICE int binarySearchVector(double a_x, Vector< double > const &a_Xs, bool a_boundIndex=false)
Definition MCGIDI.hpp:318

◆ evaluateFormFactor()

LUPI_HOST_DEVICE double MCGIDI::Distributions::CoherentPhotoAtomicScattering::evaluateFormFactor ( double a_energyIn,
double a_mu ) const

This method evaluates the coherent photo-atomic form factor at the projectile energy a_energy and product cosine of angle a_mu.

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

Definition at line 702 of file MCGIDI_distributions.cc.

702 {
703
704 double X = a_energyIn * sqrt( 0.5 * ( 1 - a_mu ) );
705 int intLowerIndex = binarySearchVector( X, m_energies );
706
707 if( intLowerIndex < 1 ) {
708 if( intLowerIndex == 0 ) return( m_formFactor[0] );
709 if( intLowerIndex == -2 ) return( m_formFactor[0] ); // This should never happend for proper a_energyIn and a_mu.
710 return( m_formFactor.back( ) );
711 }
712
713 std::size_t lowerIndex = static_cast<std::size_t>( intLowerIndex );
714
715 return( m_formFactor[lowerIndex] * pow( X / m_energies[lowerIndex] , m_a[lowerIndex] ) );
716}

Referenced by angleBiasing(), and evaluate().

◆ sample()

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

This method samples the outgoing product data from the coherent photo-atomic scattering law. It also samples the outgoing phi uniformly between 0 and 2 pi.

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

975 {
976
977 a_input.m_energyOut1 = a_X;
978
979 int intLowerIndex = binarySearchVector( a_X, m_energies );
980
981 if( intLowerIndex < 1 ) {
982 do {
983 a_input.m_mu = 1.0 - 2.0 * a_rng( );
984 } while( ( 1.0 + a_input.m_mu * a_input.m_mu ) < 2.0 * a_rng( ) ); }
985 else {
986 std::size_t lowerIndex = static_cast<std::size_t>( intLowerIndex );
987 double _a = m_a[lowerIndex];
988 double X_i = m_energies[lowerIndex];
989 double formFactor_i = m_formFactor[lowerIndex];
990 double formFactor_X_i = formFactor_i * X_i;
991 double Z = a_X / X_i;
992 double realAnomalousFactor = 0.0;
993 double imaginaryAnomalousFactor = 0.0;
994
995 if( m_anomalousDataPresent ) {
996 realAnomalousFactor = m_realAnomalousFactor->evaluate( a_X );
997 imaginaryAnomalousFactor = m_imaginaryAnomalousFactor->evaluate( a_X );
998 }
999
1000 double anomalousFactorSquared = realAnomalousFactor * realAnomalousFactor + imaginaryAnomalousFactor * imaginaryAnomalousFactor;
1001 double normalization = m_integratedFormFactorSquared[lowerIndex] + formFactor_X_i * formFactor_X_i * Z_a( Z, 2.0 * _a + 2.0 );
1002
1003 anomalousFactorSquared = 0.0;
1004 if( anomalousFactorSquared != 0.0 ) {
1005 double integratedFormFactor_i = m_integratedFormFactor[lowerIndex] + formFactor_X_i * Z_a( Z, _a + 2.0 );
1006
1007 normalization += 2.0 * integratedFormFactor_i * realAnomalousFactor + 0.5 * anomalousFactorSquared * a_X * a_X;
1008 }
1009
1010 do {
1011 double partialIntegral = a_rng( ) * normalization;
1012 double X;
1013 if( anomalousFactorSquared == 0.0 ) {
1014 intLowerIndex = binarySearchVector( partialIntegral, m_integratedFormFactorSquared );
1015 lowerIndex = static_cast<std::size_t>( intLowerIndex );
1016
1017 if( lowerIndex == 0 ) {
1018 X = sqrt( 2.0 * partialIntegral ) / m_formFactor[0]; }
1019 else {
1020 double remainer = partialIntegral - m_integratedFormFactorSquared[lowerIndex];
1021 double epsilon = 2.0 * m_a[lowerIndex] + 2.0;
1022
1023 X_i = m_energies[lowerIndex];
1024 formFactor_i = m_formFactor[lowerIndex];
1025 formFactor_X_i = formFactor_i * X_i;
1026
1027 remainer /= formFactor_X_i * formFactor_X_i;
1028 if( fabs( epsilon ) < 1e-6 ) {
1029 X = X_i * exp( remainer ); }
1030 else {
1031 X = X_i * pow( 1.0 + epsilon * remainer, 1.0 / epsilon );
1032 }
1033 } }
1034 else { // Currently not implemented.
1035 X = 0.5 * a_X;
1036 }
1037 double X_E = X / a_X;
1038 a_input.m_mu = 1.0 - 2.0 * X_E * X_E;
1039 } while( ( 1.0 + a_input.m_mu * a_input.m_mu ) < 2.0 * a_rng( ) );
1040 }
1041
1042 a_input.m_phi = 2.0 * M_PI * a_rng( );
1043 a_input.m_frame = productFrame( );
1044}
G4double epsilon(G4double density, G4double temperature)
LUPI_HOST_DEVICE GIDI::Frame productFrame() const

◆ serialize()

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

770 {
771
772 Distribution::serialize( a_buffer, a_mode );
773
774 DATA_MEMBER_INT( m_anomalousDataPresent, a_buffer, a_mode );
775 DATA_MEMBER_VECTOR_DOUBLE( m_energies, a_buffer, a_mode );
776 DATA_MEMBER_VECTOR_DOUBLE( m_formFactor, a_buffer, a_mode );
777 DATA_MEMBER_VECTOR_DOUBLE( m_a, a_buffer, a_mode );
778 DATA_MEMBER_VECTOR_DOUBLE( m_integratedFormFactor, a_buffer, a_mode );
779 DATA_MEMBER_VECTOR_DOUBLE( m_integratedFormFactorSquared, a_buffer, a_mode );
780 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm1_1, a_buffer, a_mode );
781 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm1_3, a_buffer, a_mode );
782 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm1_5, a_buffer, a_mode );
783 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm2_1, a_buffer, a_mode );
784 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm2_3, a_buffer, a_mode );
785 DATA_MEMBER_VECTOR_DOUBLE( m_probabilityNorm2_5, a_buffer, a_mode );
786
787 if( m_anomalousDataPresent ) {
788 m_realAnomalousFactor = serializeFunction1d_d1( a_buffer, a_mode, m_realAnomalousFactor );
789 m_imaginaryAnomalousFactor = serializeFunction1d_d1( a_buffer, a_mode, m_imaginaryAnomalousFactor );
790 }
791}
#define DATA_MEMBER_VECTOR_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 Functions::Function1d_d1 * serializeFunction1d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d)

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