Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
GIDI::Reaction Class Reference

#include <GIDI.hpp>

Inheritance diagram for GIDI::Reaction:

Public Member Functions

 Reaction (int a_ENDF_MT, std::string const &a_fissionGenre)
 Reaction (Construction::Settings const &a_construction, HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pops, PoPI::Database const &a_internalPoPs, Protare const &a_protare, Styles::Suite const *a_styles)
 ~Reaction ()
bool active () const
void setActive (bool a_active)
std::size_t reactionIndex () const
int depth () const
int ENDF_MT () const
int ENDL_C () const
int ENDL_S () const
std::string const & fissionGenre () const
bool isPairProduction () const
bool isPhotoAtomicIncoherentScattering () const
bool RutherfordScatteringPresent () const
bool onlyRutherfordScatteringPresent () const
bool nuclearPlusInterferencePresent () const
ComponentdoubleDifferentialCrossSection ()
Component const & doubleDifferentialCrossSection () const
ComponentcrossSection ()
Component const & crossSection () const
ComponentavailableEnergy ()
Component const & availableEnergy () const
ComponentavailableMomentum ()
Component const & availableMomentum () const
OutputChanneloutputChannel () const
void setOutputChannel (OutputChannel *a_outputChannel)
void modifiedMultiGroupElasticForTNSL (std::map< std::string, std::size_t > const &a_maximumTNSL_MultiGroupIndex)
GUPI::AncestryfindInAncestry3 (std::string const &a_item)
GUPI::Ancestry const * findInAncestry3 (std::string const &a_item) const
std::string xlinkItemKey () const
bool hasFission () const
void productIDs (std::set< std::string > &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly) const
int productMultiplicity (std::string const &a_productID) const
int maximumLegendreOrder (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
double threshold () const
double crossSectionThreshold () const
double twoBodyThreshold () const
bool areAllProductsTracked (Transporting::Particles const &a_particles) const
Vector multiGroupCrossSection (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_label="") const
Vector multiGroupQ (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, bool a_final) const
Vector multiGroupMultiplicity (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Matrix multiGroupProductMatrix (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::string const &a_productID, std::size_t a_order) const
Matrix multiGroupFissionMatrix (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::size_t a_order) const
Vector multiGroupAvailableEnergy (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
Vector multiGroupAverageEnergy (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Vector multiGroupDepositionEnergy (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles) const
Vector multiGroupAvailableMomentum (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
Vector multiGroupAverageMomentum (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Vector multiGroupDepositionMomentum (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles) const
Vector multiGroupGain (LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID, std::string const &a_projectileID) const
void delayedNeutronProducts (DelayedNeutronProducts &a_delayedNeutronProducts) const
void incompleteParticles (Transporting::Settings const &a_settings, std::set< std::string > &a_incompleteParticles) const
void continuousEnergyProductData (Transporting::Settings const &a_settings, std::string const &a_particleID, double a_energy, double &a_productEnergy, double &a_productMomentum, double &a_productGain, bool a_ignoreIncompleteParticles) const
void mapContinuousEnergyProductData (Transporting::Settings const &a_settings, std::string const &a_particleID, std::vector< double > const &a_energies, std::size_t a_offset, std::vector< double > &a_productEnergies, std::vector< double > &a_productMomenta, std::vector< double > &a_productGains, bool a_ignoreIncompleteParticles) const
bool modifyCrossSection (Functions::XYs1d const *a_offset, Functions::XYs1d const *a_slope, bool a_updateMultiGroup=false)
bool modifiedCrossSection (Functions::XYs1d const *a_offset, Functions::XYs1d const *a_slope)
void recalculateMultiGroupData (ProtareSingle const *a_protare, Styles::TemperatureInfo const &a_temperatureInfo)
void calculateMultiGroupData (ProtareSingle const *a_protare, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_heatedMultiGroupLabel, MultiGroupCalulationInformation const &a_multiGroupCalulationInformation)
void toXMLList (GUPI::WriteInfo &a_writeInfo, std::string const &a_indent="") const
Public Member Functions inherited from GIDI::Form
 Form (FormType a_type)
 Form (std::string const &a_moniker, FormType a_type, std::string const &a_label)
 Form (HAPI::Node const &a_node, SetupInfo &a_setupInfo, FormType a_type, Suite *a_suite=nullptr)
 Form (Form const &a_form)
virtual ~Form ()
Formoperator= (Form const &a_rhs)
Suiteparent () const
std::string const & label () const
void setLabel (std::string const &a_label)
virtual std::string actualMoniker () const
std::string const & keyName () const
void setKeyName (std::string const &a_keyName)
std::string const & keyValue () const
virtual void setKeyValue (std::string const &a_keyName) const
FormType type () const
Form const * sibling (std::string a_label) const
GUPI::AncestryfindInAncestry3 (LUPI_maybeUnused std::string const &a_item)
GUPI::Ancestry const * findInAncestry3 (LUPI_maybeUnused std::string const &a_item) const
Public Member Functions inherited from GUPI::Ancestry
 Ancestry (std::string const &a_moniker, std::string const &a_attribute="")
virtual ~Ancestry ()
Ancestryoperator= (Ancestry const &a_ancestry)
std::string const & moniker () const
void setMoniker (std::string const &a_moniker)
Ancestryancestor ()
Ancestry const * ancestor () const
void setAncestor (Ancestry *a_ancestor)
std::string attribute () const
Ancestryroot ()
Ancestry const * root () const
bool isChild (Ancestry *a_instance)
bool isParent (Ancestry *a_parent)
bool isRoot () const
AncestryfindInAncestry (std::string const &a_href)
Ancestry const * findInAncestry (std::string const &a_href) const
virtual LUPI_HOST void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
std::string toXLink () const
void printXML () const

Friends

class ProtareSingle

Additional Inherited Members

Static Public Member Functions inherited from GUPI::Ancestry
static std::string buildXLinkItemKey (std::string const &a_name, std::string const &a_key)

Detailed Description

Definition at line 4245 of file GIDI.hpp.

Constructor & Destructor Documentation

◆ Reaction() [1/2]

GIDI::Reaction::Reaction ( int a_ENDF_MT,
std::string const & a_fissionGenre )

Parses a <reaction> node.

Definition at line 21 of file GIDI_reaction.cc.

21 :
23 m_reactionIndex( 0 ),
24 m_active( true ),
25 m_ENDF_MT( a_ENDF_MT ),
26 m_fissionGenre( a_fissionGenre ),
27 m_QThreshold( 0.0 ),
28 m_crossSectionThreshold( 0.0 ),
29 m_twoBodyThreshold( 0.0 ),
30 m_isPairProduction( false ),
31 m_isPhotoAtomicIncoherentScattering( false ),
32 m_RutherfordScatteringPresent( false ),
33 m_onlyRutherfordScatteringPresent( false ),
34 m_nuclearPlusInterferencePresent( false ),
35 m_decayPositronium( false ),
36 m_doubleDifferentialCrossSection( GIDI_doubleDifferentialCrossSectionChars, GIDI_labelChars ),
37 m_crossSection( GIDI_crossSectionChars, GIDI_labelChars ),
38 m_availableEnergy( GIDI_availableEnergyChars, GIDI_labelChars ),
39 m_availableMomentum( GIDI_availableMomentumChars, GIDI_labelChars ),
40 m_outputChannel( ) {
41
43
44 ENDL_CFromENDF_MT( m_ENDF_MT, &m_ENDL_C, &m_ENDL_S );
45}
#define GIDI_crossSectionChars
Definition GIDI.hpp:215
#define GIDI_availableMomentumChars
Definition GIDI.hpp:217
#define GIDI_reactionChars
Definition GIDI.hpp:180
#define GIDI_availableEnergyChars
Definition GIDI.hpp:216
#define GIDI_labelChars
Definition GIDI.hpp:438
#define GIDI_doubleDifferentialCrossSectionChars
Definition GIDI.hpp:214
Form(FormType a_type)
Definition GIDI_form.cc:25
void setMoniker(std::string const &a_moniker)
Definition GUPI.hpp:103
int ENDL_CFromENDF_MT(int ENDF_MT, int *ENDL_C, int *ENDL_S)

◆ Reaction() [2/2]

GIDI::Reaction::Reaction ( Construction::Settings const & a_construction,
HAPI::Node const & a_node,
SetupInfo & a_setupInfo,
PoPI::Database const & a_pops,
PoPI::Database const & a_internalPoPs,
Protare const & a_protare,
Styles::Suite const * a_styles )

Parses a <reaction> node.

Parameters
a_construction[in] Used to pass user options to the constructor.
a_node[in] The reaction HAPI::Node to be parsed and used to construct the reaction.
a_setupInfo[in] Information create my the Protare constructor to help in parsing.
a_pops[in] The external PoPI::Database instance used to get particle indices and possibly other particle information.
a_internalPoPs[in] The internal PoPI::Database instance used to get particle indices and possibly other particle information. This is the <PoPs> node under the <reactionSuite> node.
a_protare[in] The GIDI::Protare this reaction belongs to.
a_styles[in] The <styles> node under the <reactionSuite> node.

Definition at line 60 of file GIDI_reaction.cc.

61 :
62 Form( a_node, a_setupInfo, FormType::reaction ),
63 m_active( true ),
64 m_ENDF_MT( a_node.attribute_as_int( GIDI_ENDF_MT_Chars ) ),
65 m_ENDL_C( 0 ),
66 m_ENDL_S( 0 ),
67 m_fissionGenre( a_node.attribute_as_string( GIDI_fissionGenreChars ) ),
68 m_QThreshold( 0.0 ),
69 m_crossSectionThreshold( 0.0 ),
70 m_twoBodyThreshold( 0.0 ),
71 m_isPairProduction( false ),
72 m_isPhotoAtomicIncoherentScattering( false ),
73 m_RutherfordScatteringPresent( false ),
74 m_onlyRutherfordScatteringPresent( false ),
75 m_nuclearPlusInterferencePresent( false ),
76 m_decayPositronium( false ),
77 m_doubleDifferentialCrossSection( a_construction, GIDI_doubleDifferentialCrossSectionChars, GIDI_labelChars, a_node, a_setupInfo, a_pops,
78 a_internalPoPs, parseDoubleDifferentialCrossSectionSuite, a_styles ),
79 m_crossSection( a_construction, GIDI_crossSectionChars, GIDI_labelChars, a_node, a_setupInfo, a_pops, a_internalPoPs, parseCrossSectionSuite, a_styles ),
80 m_availableEnergy( a_construction, GIDI_availableEnergyChars, GIDI_labelChars, a_node, a_setupInfo, a_pops, a_internalPoPs, parseAvailableSuite, a_styles ),
81 m_availableMomentum( a_construction, GIDI_availableMomentumChars, GIDI_labelChars, a_node, a_setupInfo, a_pops, a_internalPoPs, parseAvailableSuite, a_styles ),
82 m_outputChannel( nullptr ) {
83
84 m_isPairProduction = label( ).find( "pair production" ) != std::string::npos;
85 m_decayPositronium = m_isPairProduction && a_construction.decayPositronium( );
86 m_isPhotoAtomicIncoherentScattering = false;
87 if( m_doubleDifferentialCrossSection.size( ) > 0 )
88 m_isPhotoAtomicIncoherentScattering = m_doubleDifferentialCrossSection.get<Form>( 0 )->type( ) == FormType::incoherentPhotonScattering;
89
90 m_doubleDifferentialCrossSection.setAncestor( this );
91 m_crossSection.setAncestor( this );
92 m_availableEnergy.setAncestor( this );
93 m_availableMomentum.setAncestor( this );
94
95 a_setupInfo.m_outputChannelLevel = 0;
96 m_outputChannel = new OutputChannel( a_construction, a_node.child( GIDI_outputChannelChars ), a_setupInfo, a_pops,
97 a_internalPoPs, a_styles, hasFission( ), true );
98 m_outputChannel->setAncestor( this );
99
100 HAPI::Node const CoulombPlusNuclearElastic = a_node.child( GIDI_doubleDifferentialCrossSectionChars ).first_child( );
101 if( CoulombPlusNuclearElastic.name( ) == GIDI_CoulombPlusNuclearElasticChars ) { // Check RutherfordScattering.
102 m_RutherfordScatteringPresent = true;
103 m_onlyRutherfordScatteringPresent = true;
104 for( HAPI::Node child = CoulombPlusNuclearElastic.first_child( ); !child.empty( ); child.to_next_sibling( ) ) {
105 if( child.name( ) != GIDI_RutherfordScatteringChars ) m_onlyRutherfordScatteringPresent = false;
106 if( child.name( ) == GIDI_nuclearPlusInterferenceChars ) m_nuclearPlusInterferencePresent = true;
107 }
108 }
109
110 ENDL_CFromENDF_MT( m_ENDF_MT, &m_ENDL_C, &m_ENDL_S );
111 if( a_setupInfo.m_isENDL_C_9 ) m_ENDL_C = 9;
112 if( m_ENDL_C == 20 ) {
113 Product *product = m_outputChannel->products( ).get<Product>( 0 );
114 if( product->particle( ).ID( ) == "H1" ) m_ENDL_C = 21;
115 }
116 else if( m_ENDL_C == 22 ) {
117 Product *product = m_outputChannel->products( ).get<Product>( 0 );
118 if( product->particle( ).ID( ) == "H2" ) m_ENDL_C = 35;
119 }
120
121 if( ( a_construction.parseMode( ) != Construction::ParseMode::outline ) && ( a_construction.parseMode( ) != Construction::ParseMode::readOnly ) ) {
122 double _Q = 0.0;
123 Form const &QForm = *(m_outputChannel->Q( ).get<Form>( 0 ));
124 switch( QForm.type( ) ) {
125 case FormType::constant1d : {
126 Functions::Constant1d const &constant1d = static_cast<Functions::Constant1d const &>( QForm );
127 _Q = constant1d.value( );
128 break; }
129 case FormType::XYs1d : {
130 Functions::XYs1d const &xys1d = static_cast<Functions::XYs1d const &>( QForm );
131 _Q = xys1d.evaluate( 0.0 );
132 break; }
133 case FormType::gridded1d : { // Should be a special reaction (e.g., summed multi-group data) and can be ignored.
134 _Q = 0.0;
135 break; }
136 default :
137 throw Exception( "Reaction::Reaction: unsupported Q form " + QForm.label( ) );
138 }
139 m_twoBodyThreshold = -_Q * a_protare.thresholdFactor( );
140 _Q *= -1;
141 if( _Q <= 0.0 ) _Q = 0.0;
142 m_QThreshold = a_protare.thresholdFactor( ) * _Q;
143
144 if( _Q > 0.0 ) { // Try to do a better job determining m_crossSectionThreshold.
145 std::vector<Suite::const_iterator> monikers = a_protare.styles( ).findAllOfMoniker( GIDI_griddedCrossSectionStyleChars );
146 if( ( a_construction.parseMode( ) != Construction::ParseMode::multiGroupOnly ) && ( monikers.size( ) > 0 ) ) {
147 Styles::GriddedCrossSection const &griddedCrossSection = static_cast<Styles::GriddedCrossSection const &>( **monikers[0] );
148 Grid grid = griddedCrossSection.grid( );
149
150 if( m_crossSection.has( griddedCrossSection.label( ) ) ) {
151 Functions::Ys1d const &ys1d = static_cast<Functions::Ys1d const &>( *m_crossSection.get<Functions::Ys1d>( griddedCrossSection.label( ) ) );
152 m_crossSectionThreshold = grid[ys1d.start( )];
153 }
154 }
155 if( m_crossSectionThreshold == 0.0 ) { // Should also check 'evaluate' style before using m_QThreshold as a default.
156 m_crossSectionThreshold = m_QThreshold;
157 }
158 }
159 if( m_twoBodyThreshold > 0.0 ) m_twoBodyThreshold = m_crossSectionThreshold;
160 }
161}
#define GIDI_nuclearPlusInterferenceChars
Definition GIDI.hpp:377
#define GIDI_ENDF_MT_Chars
Definition GIDI.hpp:420
#define GIDI_griddedCrossSectionStyleChars
Definition GIDI.hpp:251
#define GIDI_outputChannelChars
Definition GIDI.hpp:227
#define GIDI_RutherfordScatteringChars
Definition GIDI.hpp:376
#define GIDI_fissionGenreChars
Definition GIDI.hpp:410
#define GIDI_CoulombPlusNuclearElasticChars
Definition GIDI.hpp:375
std::string const & label() const
Definition GIDI.hpp:658
FormType type() const
Definition GIDI.hpp:667
bool hasFission() const
std::string name() const
Definition HAPI_Node.cc:136
bool empty() const
Definition HAPI_Node.cc:150
Node first_child() const
Definition HAPI_Node.cc:82
Form * parseCrossSectionSuite(Construction::Settings const &a_construction, Suite *parent, HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pop, PoPI::Database const &a_internalPoPs, std::string const &a_name, Styles::Suite const *a_styles)
@ incoherentPhotonScattering
Definition GIDI.hpp:135
Form * parseAvailableSuite(Construction::Settings const &a_construction, Suite *a_parent, HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pops, PoPI::Database const &a_internalPoPs, std::string const &a_name, Styles::Suite const *a_styles)
Form * parseDoubleDifferentialCrossSectionSuite(Construction::Settings const &a_construction, Suite *a_parent, HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pops, PoPI::Database const &a_internalPoPs, std::string const &a_name, Styles::Suite const *a_styles)

◆ ~Reaction()

GIDI::Reaction::~Reaction ( )

Definition at line 166 of file GIDI_reaction.cc.

166 {
167
168 if( m_outputChannel != nullptr ) delete m_outputChannel;
169}

Member Function Documentation

◆ active()

◆ areAllProductsTracked()

bool GIDI::Reaction::areAllProductsTracked ( Transporting::Particles const & a_particles) const

Returns true if all outgoing particles (i.e., products) are specifed in a_particles. That is, the user will be tracking all products of this reaction.

Parameters
a_particles[in] The list of particles to be transported.
Returns
bool.

Definition at line 331 of file GIDI_reaction.cc.

331 {
332
333 if( hasFission( ) || m_isPhotoAtomicIncoherentScattering ) return( false );
334
335 return( m_outputChannel->areAllProductsTracked( a_particles ) );
336}

Referenced by GIDI::ProtareSingle::multiGroupDepositionEnergy(), and multiGroupDepositionEnergy().

◆ availableEnergy() [1/2]

Component & GIDI::Reaction::availableEnergy ( )
inline

Returns a reference to the m_availableEnergy member.

Definition at line 4302 of file GIDI.hpp.

Referenced by multiGroupDepositionEnergy().

◆ availableEnergy() [2/2]

Component const & GIDI::Reaction::availableEnergy ( ) const
inline

Returns a reference to the m_availableEnergy member.

Definition at line 4303 of file GIDI.hpp.

◆ availableMomentum() [1/2]

Component & GIDI::Reaction::availableMomentum ( )
inline

Returns a reference to the m_availableMomentum member.

Definition at line 4304 of file GIDI.hpp.

◆ availableMomentum() [2/2]

Component const & GIDI::Reaction::availableMomentum ( ) const
inline

Returns a reference to the m_availableMomentum member.

Definition at line 4305 of file GIDI.hpp.

◆ calculateMultiGroupData()

void GIDI::Reaction::calculateMultiGroupData ( ProtareSingle const * a_protare,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_heatedMultiGroupLabel,
MultiGroupCalulationInformation const & a_multiGroupCalulationInformation )

This methods calculates multi-group data for all needed components and adds each component's multi-group with label a_heatedMultiGroupLabel.

Parameters
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_heatedMultiGroupLabel[in] The label of the style for the multi-group data being added.
a_multiGroupCalulationInformation[in] Store multi-group boundary and flux data used for multi-grouping.
a_crossSectionXYs1d[in[ The cross section weight.

Definition at line 867 of file GIDI_reaction.cc.

868 {
869
870 Styles::Suite const &styles = a_protare->styles( );
871 Transporting::MultiGroup multiGroup = a_multiGroupCalulationInformation.m_multiGroup;
872
873 std::vector<double> flatValues { multiGroup[0], 1.0, multiGroup[multiGroup.size( ) - 1], 1.0 };
874 Functions::XYs1d flatFunction( Axes( ), ptwXY_interpolationLinLin, flatValues );
875
876 Functions::XYs1d const *crossSectionXYs1d = static_cast<Functions::XYs1d *>(
877 m_crossSection.findInstanceOfTypeInLineage( styles, a_temperatureInfo.heatedCrossSection( ), GIDI_XYs1dChars ) );
878
879 calculate1dMultiGroupDataInComponent( a_protare, a_heatedMultiGroupLabel, a_multiGroupCalulationInformation, m_crossSection, flatFunction );
880 calculate1dMultiGroupDataInComponent( a_protare, a_heatedMultiGroupLabel, a_multiGroupCalulationInformation, m_availableEnergy, *crossSectionXYs1d );
881 calculate1dMultiGroupDataInComponent( a_protare, a_heatedMultiGroupLabel, a_multiGroupCalulationInformation, m_availableMomentum, *crossSectionXYs1d );
882
883 if( m_outputChannel != nullptr )
884 m_outputChannel->calculateMultiGroupData( a_protare, a_temperatureInfo, a_heatedMultiGroupLabel, a_multiGroupCalulationInformation, *crossSectionXYs1d );
885}
#define GIDI_XYs1dChars
Definition GIDI.hpp:288
void calculate1dMultiGroupDataInComponent(ProtareSingle const *a_protare, std::string const &a_heatedMultiGroupLabel, MultiGroupCalulationInformation const &a_multiGroupCalulationInformation, Component &a_component, Functions::XYs1d const &a_crossSection)
@ ptwXY_interpolationLinLin
Definition ptwXY.h:37

Referenced by recalculateMultiGroupData().

◆ continuousEnergyProductData()

void GIDI::Reaction::continuousEnergyProductData ( Transporting::Settings const & a_settings,
std::string const & a_particleID,
double a_energy,
double & a_productEnergy,
double & a_productMomentum,
double & a_productGain,
bool a_ignoreIncompleteParticles ) const

Returns, via arguments, the average energy and momentum, and gain for product with particle id a_particleID.

Parameters
a_settings[in] Specifies the requested label.
a_particleID[in] The particle id of the product.
a_energy[in] The energy of the projectile.
a_productEnergy[in] The average energy of the product.
a_productMomentum[in] The average momentum of the product.
a_productGain[in] The gain of the product.
a_ignoreIncompleteParticles[in] If true, incomplete particles are ignore, otherwise a throw is executed.

Definition at line 666 of file GIDI_reaction.cc.

667 {
668
669 a_productEnergy = 0.0;
670 a_productMomentum = 0.0;
671 a_productGain = 0.0;
672
673// if( ENDF_MT( ) == 516 ) return; // FIXME, may be something wrong with the way FUDGE converts ENDF to GNDS.
674
675 if( m_outputChannel != nullptr )
676 m_outputChannel->continuousEnergyProductData( a_settings, a_particleID, a_energy, a_productEnergy, a_productMomentum, a_productGain,
677 a_ignoreIncompleteParticles );
678}

◆ crossSection() [1/2]

Component & GIDI::Reaction::crossSection ( )
inline

Returns a reference to the m_crossSection member.

Definition at line 4299 of file GIDI.hpp.

◆ crossSection() [2/2]

Component const & GIDI::Reaction::crossSection ( ) const
inline

Returns a reference to the m_crossSection member.

Definition at line 4300 of file GIDI.hpp.

◆ crossSectionThreshold()

double GIDI::Reaction::crossSectionThreshold ( ) const
inline

Returns the value of the m_crossSectionThreshold member.

Definition at line 4324 of file GIDI.hpp.

◆ delayedNeutronProducts()

void GIDI::Reaction::delayedNeutronProducts ( DelayedNeutronProducts & a_delayedNeutronProducts) const

Appends a DelayedNeutronProduct instance for each delayed neutron in m_delayedNeutrons.

Parameters
a_delayedNeutronProducts[in/out] The list to append the delayed neutrons to.

Definition at line 636 of file GIDI_reaction.cc.

636 {
637
638 if( m_outputChannel != nullptr ) m_outputChannel->delayedNeutronProducts( a_delayedNeutronProducts );
639}

Referenced by GIDI::ProtareSingle::delayedNeutronProducts().

◆ depth()

int GIDI::Reaction::depth ( ) const
inline

Returns the maximum product depth for this reaction.

Definition at line 4283 of file GIDI.hpp.

◆ doubleDifferentialCrossSection() [1/2]

Component & GIDI::Reaction::doubleDifferentialCrossSection ( )
inline

Returns a reference to the m_doubleDifferentialCrossSection member.

Definition at line 4297 of file GIDI.hpp.

Referenced by MCGIDI::Distributions::parseGIDI().

◆ doubleDifferentialCrossSection() [2/2]

Component const & GIDI::Reaction::doubleDifferentialCrossSection ( ) const
inline

Returns a reference to the m_doubleDifferentialCrossSection member.

Definition at line 4298 of file GIDI.hpp.

◆ ENDF_MT()

int GIDI::Reaction::ENDF_MT ( ) const
inline

Returns the value of the m_ENDF_MT member.

Definition at line 4284 of file GIDI.hpp.

Referenced by toXMLList().

◆ ENDL_C()

int GIDI::Reaction::ENDL_C ( ) const
inline

Returns the value of the m_ENDL_C member.

Definition at line 4285 of file GIDI.hpp.

Referenced by GIDI::Protare::reactionIndicesMatchingENDLCValues().

◆ ENDL_S()

int GIDI::Reaction::ENDL_S ( ) const
inline

Returns the value of the m_ENDL_S member.

Definition at line 4286 of file GIDI.hpp.

◆ findInAncestry3() [1/2]

GUPI::Ancestry * GIDI::Reaction::findInAncestry3 ( std::string const & a_item)
virtual

Used by GUPI::Ancestry to tranverse GNDS nodes. This method returns a pointer to a derived class' a_item member or nullptr if none exists.

Parameters
a_item[in] The name of the class member whose pointer is to be return.
Returns
The pointer to the class member or nullptr if class does not have a member named a_item.

Implements GUPI::Ancestry.

Definition at line 292 of file GIDI_reaction.cc.

292 {
293
294 if( a_item == GIDI_doubleDifferentialCrossSectionChars ) return( &m_doubleDifferentialCrossSection );
295 if( a_item == GIDI_crossSectionChars ) return( &m_crossSection );
296 if( a_item == GIDI_availableEnergyChars ) return( &m_availableEnergy );
297 if( a_item == GIDI_availableMomentumChars ) return( &m_availableMomentum );
298 if( a_item == GIDI_outputChannelChars ) return( m_outputChannel );
299
300 return( nullptr );
301}

◆ findInAncestry3() [2/2]

GUPI::Ancestry const * GIDI::Reaction::findInAncestry3 ( std::string const & a_item) const
virtual

Used by GUPI::Ancestry to tranverse GNDS nodes. This method returns a pointer to a derived class' a_item member or nullptr if none exists.

Parameters
a_item[in] The name of the class member whose pointer is to be return.
Returns
The pointer to the class member or nullptr if class does not have a member named a_item.

Implements GUPI::Ancestry.

Definition at line 311 of file GIDI_reaction.cc.

311 {
312
313 if( a_item == GIDI_doubleDifferentialCrossSectionChars ) return( &m_doubleDifferentialCrossSection );
314 if( a_item == GIDI_crossSectionChars ) return( &m_crossSection );
315 if( a_item == GIDI_availableEnergyChars ) return( &m_availableEnergy );
316 if( a_item == GIDI_availableMomentumChars ) return( &m_availableMomentum );
317 if( a_item == GIDI_outputChannelChars ) return( m_outputChannel );
318
319 return( nullptr );
320}

◆ fissionGenre()

std::string const & GIDI::Reaction::fissionGenre ( ) const
inline

Definition at line 4287 of file GIDI.hpp.

4287{ return( m_fissionGenre ); }

◆ hasFission()

bool GIDI::Reaction::hasFission ( ) const

Returns true if at least one output channel contains a fission channel.

Returns
true if at least one output channel is a fission channel.

Definition at line 246 of file GIDI_reaction.cc.

246 {
247
248 if( m_ENDF_MT == 18 ) return( true );
249 if( m_ENDF_MT == 19 ) return( true );
250 if( m_ENDF_MT == 20 ) return( true );
251 if( m_ENDF_MT == 21 ) return( true );
252 if( m_ENDF_MT == 38 ) return( true );
253
254 return( false );
255}

Referenced by areAllProductsTracked(), GIDI::ProtareSingle::delayedNeutronProducts(), GIDI::ProtareSingle::hasFission(), GIDI::ProtareSingle::isDelayedFissionNeutronComplete(), GIDI::ProtareSingle::multiGroupFissionGammaMultiplicity(), GIDI::ProtareSingle::multiGroupFissionMatrix(), multiGroupFissionMatrix(), GIDI::ProtareSingle::multiGroupFissionNeutronMultiplicity(), and Reaction().

◆ incompleteParticles()

void GIDI::Reaction::incompleteParticles ( Transporting::Settings const & a_settings,
std::set< std::string > & a_incompleteParticles ) const

If this has an output channel, this output channel's **incompleteParticles* is called.

Parameters
a_settings[in] Specifies the requested label.
a_incompleteParticles[out] The list of particles whose completeParticle method returns false.

Definition at line 648 of file GIDI_reaction.cc.

648 {
649
650 if( m_outputChannel != nullptr ) m_outputChannel->incompleteParticles( a_settings, a_incompleteParticles );
651 if( m_isPairProduction ) a_incompleteParticles.erase( PoPI::IDs::photon ); // Kludge for old processed files.
652}
static std::string const photon
Definition PoPI.hpp:162

Referenced by GIDI::ProtareSingle::incompleteParticles().

◆ isPairProduction()

bool GIDI::Reaction::isPairProduction ( ) const
inline

Returns the value of the m_isPairProduction member.

Definition at line 4288 of file GIDI.hpp.

◆ isPhotoAtomicIncoherentScattering()

bool GIDI::Reaction::isPhotoAtomicIncoherentScattering ( ) const
inline

Returns the value of the m_isPhotoAtomicIncoherentScattering member.

Definition at line 4289 of file GIDI.hpp.

◆ mapContinuousEnergyProductData()

void GIDI::Reaction::mapContinuousEnergyProductData ( Transporting::Settings const & a_settings,
std::string const & a_particleID,
std::vector< double > const & a_energies,
std::size_t a_offset,
std::vector< double > & a_productEnergies,
std::vector< double > & a_productMomenta,
std::vector< double > & a_productGains,
bool a_ignoreIncompleteParticles ) const

Modifies the average product energies, momenta and gains for product with particle id a_particleID.

Parameters
a_settings[in] Specifies user options.
a_particleID[in] The particle id of the product.
a_energies[in] The vector of energies to map the data to.
a_offset[in] The index of the first energy whose data are to be added to the vectors.
a_productEnergies[out] The vector of average energies of the product.
a_productMomenta[out] The vector of average momenta of the product.
a_productGains[out] The vector of gain of the product.
a_ignoreIncompleteParticles[in] If true, incomplete particles are ignore, otherwise a throw is executed.

Definition at line 693 of file GIDI_reaction.cc.

695 {
696
697// if( ENDF_MT( ) == 516 ) return; // FIXME, may be something wrong with the way FUDGE converts ENDF to GNDS.
698
699 if( m_outputChannel != nullptr )
700 m_outputChannel->mapContinuousEnergyProductData( a_settings, a_particleID, a_energies, a_offset, a_productEnergies,
701 a_productMomenta, a_productGains, a_ignoreIncompleteParticles );
702}

◆ maximumLegendreOrder()

int GIDI::Reaction::maximumLegendreOrder ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_productID ) const

Determines the maximum Legredre order present in the multi-group transfer matrix for a give product for a give label. Inspects all products produced by this reaction.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_productID[in] Particle id of the requested product.
Returns
The maximum Legredre order. If no transfer matrix data are present for the requested product, -1 is returned.

Definition at line 208 of file GIDI_reaction.cc.

209 {
210
211 if( m_isPairProduction && ( a_productID == PoPI::IDs::photon ) ) return( 0 );
212 return( m_outputChannel->maximumLegendreOrder( a_smr, a_settings, a_temperatureInfo, a_productID ) );
213}

Referenced by GIDI::ProtareSingle::maximumLegendreOrder().

◆ modifiedCrossSection()

bool GIDI::Reaction::modifiedCrossSection ( Functions::XYs1d const * a_offset,
Functions::XYs1d const * a_slope )

Thid methid is deprecated, use modifyCrossSection instead. See modifyCrossSection for useage.

Definition at line 820 of file GIDI_reaction.cc.

820 {
821
822 return modifyCrossSection( a_offset, a_slope, false );
823}
bool modifyCrossSection(Functions::XYs1d const *a_offset, Functions::XYs1d const *a_slope, bool a_updateMultiGroup=false)

◆ modifiedMultiGroupElasticForTNSL()

void GIDI::Reaction::modifiedMultiGroupElasticForTNSL ( std::map< std::string, std::size_t > const & a_maximumTNSL_MultiGroupIndex)

Only for internal use. Called by ProtareTNSL instance to zero the lower energy multi-group data covered by the TNSL ProtareSingle.

Parameters
a_maximumTNSL_MultiGroupIndex[in] A map that contains labels for heated multi-group data and the last valid group boundary for the TNSL data for that boundary.

Definition at line 276 of file GIDI_reaction.cc.

276 {
277
278 m_crossSection.modifiedMultiGroupElasticForTNSL( a_maximumTNSL_MultiGroupIndex );
279 m_availableEnergy.modifiedMultiGroupElasticForTNSL( a_maximumTNSL_MultiGroupIndex );
280 m_availableMomentum.modifiedMultiGroupElasticForTNSL( a_maximumTNSL_MultiGroupIndex );
281 m_outputChannel->modifiedMultiGroupElasticForTNSL( a_maximumTNSL_MultiGroupIndex );
282}

◆ modifyCrossSection()

bool GIDI::Reaction::modifyCrossSection ( Functions::XYs1d const * a_offset,
Functions::XYs1d const * a_slope,
bool a_updateMultiGroup = false )

This method modifies the cross section for this reaction as

     crossSection = a_offset + a_slope * crossSection

Either or both of a_offset and a_slope can be an empty Functions::XYs1d instance or a nullptr. If both a_offset and a_slope are non-empty Functions::XYs1d instances, the domains of both must be the same. If the returned value is false, no data are changed.

Parameters
a_offset[in] A pointer to a XYs1d function for the offset.
a_slope[in] A pointer to a XYs1d function for the slope.
a_updateMultiGroup[in] If true, the multi-group data are also modified.
Returns
true if data can be modified and false otherwise.

Definition at line 720 of file GIDI_reaction.cc.

720 {
721/*
722 -) FIXME, this does not modify the total cross section for the protare! Does it needed to be?
723*/
724
725 ProtareSingle &protare( dynamic_cast<ProtareSingle &>( *root( ) ) );
726 Styles::Suite const &styles = protare.styles( );
727 Functions::XYs1d const *offset1 = a_offset, *slope1 = a_slope;
728
729 if( a_offset == nullptr ) offset1 = new Functions::XYs1d( ); // Handle nullptr cases.
730 if( a_slope == nullptr ) slope1 = new Functions::XYs1d( );
731
732 if( offset1->size( ) == 0 ) { // Handle empty functions.
733 if( slope1->size( ) == 0 ) {
734 if( a_offset == nullptr ) delete offset1;
735 if( a_slope == nullptr ) delete slope1;
736 return( false );
737 }
738 Functions::XYs1d const *offset2 = Functions::XYs1d::makeConstantXYs1d( offset1->axes( ), slope1->domainMin( ), slope1->domainMax( ), 0.0 );
739 if( a_offset == nullptr ) delete offset1;
740 offset1 = offset2; }
741 else if( slope1->size( ) == 0 ) {
742 Functions::XYs1d const *slope2 = Functions::XYs1d::makeConstantXYs1d( slope1->axes( ), offset1->domainMin( ), offset1->domainMax( ), 1.0 );
743 if( a_slope == nullptr ) delete slope1;
744 slope1 = slope2;
745 }
746
747 double domainMin = offset1->domainMin( ), domainMax = offset1->domainMax( );
748
749 if( domainMin != slope1->domainMin( ) ) throw Exception( "GIDI::Reaction::modifyCrossSection: offset and slope domainMins differ." );
750 if( domainMax != slope1->domainMax( ) ) throw Exception( "GIDI::Reaction::modifyCrossSection: offset and slope domainMaxs differ." );
751
752 Styles::TemperatureInfos temperatureInfos = protare.temperatures( );
753
754 Functions::XYs1d *xys1d = static_cast<Functions::XYs1d *>(
755 m_crossSection.findInstanceOfTypeInLineage( styles, temperatureInfos[0].heatedCrossSection( ), GIDI_XYs1dChars ) );
756 if( xys1d == nullptr ) throw Exception( "GIDI::Reaction::modifyCrossSection: could not find XYs1d cross section." );
757
758 if( ( xys1d->domainMin( ) >= domainMax ) || ( xys1d->domainMax( ) <= domainMin ) ) {
759 if( a_offset == nullptr ) delete offset1;
760 if( a_slope == nullptr ) delete slope1;
761 return( false );
762 }
763
764 double crossSectionDomainMax = xys1d->domainMax( );
765 if( xys1d->domainMin( ) > domainMin ) domainMin = xys1d->domainMin( );
766 if( crossSectionDomainMax < domainMax ) domainMax = crossSectionDomainMax;
767
768 Functions::XYs1d offset = offset1->domainSlice( domainMin, domainMax, true );
769 Functions::XYs1d slope = slope1->domainSlice( domainMin, domainMax, true );
770 if( a_offset == nullptr ) delete offset1;
771 if( a_slope == nullptr ) delete slope1;
772
773 for( auto temperatureInfo = temperatureInfos.begin( ); temperatureInfo != temperatureInfos.end( ); ++temperatureInfo ) {
774 xys1d = static_cast<Functions::XYs1d *>(
775 m_crossSection.findInstanceOfTypeInLineage( styles, temperatureInfo->heatedCrossSection( ), GIDI_XYs1dChars ) );
776 Functions::XYs1d xys1dSliced = xys1d->domainSlice( domainMin, domainMax, true );
777 if( xys1dSliced.interpolation( ) != ptwXY_interpolationLinLin ) {
778 Functions::XYs1d *temp = xys1dSliced.asXYs1d( true, 1e-3, 1e-6, 1e-6 );
779 xys1dSliced = (*temp);
780 delete temp;
781 }
782 Functions::XYs1d modified = offset + slope * xys1dSliced;
783
784 int64_t index1;
785 ptwXYPoint *p1;
786 ptwXYPoints *ptwXY = xys1d->ptwXY( );
787 int64_t length = ptwXY_length( nullptr, ptwXY );
788 for( index1 = 0, p1 = ptwXY->points; index1 < length; ++index1, ++p1 ) {
789 if( p1->x < domainMin ) continue;
790 if( p1->x >= domainMax ) {
791 if( ( p1->x != domainMax ) || ( p1->x != crossSectionDomainMax ) ) break;
792 }
793 p1->y = modified.evaluate( p1->x );
794 }
795
796 Functions::Ys1d *ys1d = m_crossSection.get<Functions::Ys1d>( temperatureInfo->griddedCrossSection( ) );
797 std::vector<double> &Ys = ys1d->Ys( );
798 Styles::GriddedCrossSection const griddedCrossSection = *styles.get<Styles::GriddedCrossSection>( temperatureInfo->griddedCrossSection( ) );
799 nf_Buffer<double> const &grid = griddedCrossSection.grid( ).values( );
800 std::size_t start = ys1d->start( ), size = ys1d->size( ), index2;
801 for( index2 = 0; index2 < size; ++index2, ++start ) {
802 double xValue = grid[start];
803
804 if( xValue < domainMin ) continue;
805 if( xValue >= domainMax ) {
806 if( ( xValue != domainMax ) || ( xValue != crossSectionDomainMax ) ) break;
807 }
808 Ys[index2] = modified.evaluate( xValue );
809 }
810 if( a_updateMultiGroup ) recalculateMultiGroupData( &protare, *temperatureInfo );
811 }
812
813 return( true );
814}
G4ThreadLocal T * G4GeomSplitter< T >::offset
static XYs1d * makeConstantXYs1d(Axes const &a_axes, double a_domainMin, double a_domainMax, double a_value)
friend class ProtareSingle
Definition GIDI.hpp:4247
void recalculateMultiGroupData(ProtareSingle const *a_protare, Styles::TemperatureInfo const &a_temperatureInfo)
Ancestry * root()
std::vector< Styles::TemperatureInfo > TemperatureInfos
Definition GIDI.hpp:3440
struct ptwXYPoints_s ptwXYPoints
int64_t ptwXY_length(statusMessageReporting *smr, ptwXYPoints *ptwXY)
Definition ptwXY_core.c:793
struct ptwXYPoint_s ptwXYPoint
double y
Definition ptwXY.h:64
double x
Definition ptwXY.h:64
ptwXYPoint * points
Definition ptwXY.h:93

Referenced by modifiedCrossSection().

◆ multiGroupAvailableEnergy()

Vector GIDI::Reaction::multiGroupAvailableEnergy ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo ) const

Returns the multi-group, total available energy for the requested label. This is a cross section weighted available energy.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
Returns
The requested multi-group available energy as a GIDI::Vector.

Definition at line 455 of file GIDI_reaction.cc.

456 {
457
458 Vector vector( 0 );
459
460 Functions::Gridded1d const *form = dynamic_cast<Functions::Gridded1d const*>( a_settings.form( a_smr, m_availableEnergy, a_temperatureInfo, "available energy" ) );
461 if( form != nullptr ) vector = form->data( );
462
463 if( m_decayPositronium )
464 vector -= m_outputChannel->multiGroupQ( a_smr, a_settings, a_temperatureInfo, false );
465
466 return( vector );
467}

Referenced by GIDI::ProtareSingle::multiGroupAvailableEnergy(), and multiGroupDepositionEnergy().

◆ multiGroupAvailableMomentum()

Vector GIDI::Reaction::multiGroupAvailableMomentum ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo ) const

Returns the multi-group, total available momentum for the requested label. This is a cross section weighted available momentum.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
Returns
The requested multi-group available momentum as a GIDI::Vector.

Definition at line 544 of file GIDI_reaction.cc.

545 {
546
547 Vector vector( 0 );
548
549 Functions::Gridded1d const *form = dynamic_cast<Functions::Gridded1d const*>( a_settings.form( a_smr, m_availableMomentum, a_temperatureInfo, "available momentum" ) );
550 if( form != nullptr ) vector = form->data( );
551
552 return( vector );
553}

Referenced by GIDI::ProtareSingle::multiGroupAvailableMomentum(), and multiGroupDepositionMomentum().

◆ multiGroupAverageEnergy()

Vector GIDI::Reaction::multiGroupAverageEnergy ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_productID ) const

Returns the multi-group, total average energy for the requested label for the requested product. This is a cross section weighted average energy summed over all products for this reaction.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_productID[in] Particle id for the requested product.
Returns
The requested multi-group average energy as a GIDI::Vector.

Definition at line 481 of file GIDI_reaction.cc.

482 {
483
484 Vector vector( 0 );
485
486 if( m_decayPositronium ) {
487 if( a_productID == PoPI::IDs::photon ) vector += multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo ) * 2.0 * PoPI_electronMass_MeV_c2; }
488 else {
489 vector += m_outputChannel->multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, a_productID );
490 }
491
492 return( vector );
493}
#define PoPI_electronMass_MeV_c2
Definition PoPI.hpp:31
Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_label="") const

Referenced by GIDI::ProtareSingle::multiGroupAverageEnergy(), and multiGroupDepositionEnergy().

◆ multiGroupAverageMomentum()

Vector GIDI::Reaction::multiGroupAverageMomentum ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_productID ) const

Returns the multi-group, total average momentum for the requested label for the requested product. This is a cross section weighted average momentum.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_productID[in] Particle id for the requested product.
Returns
The requested multi-group average momentum as a GIDI::Vector.

Definition at line 565 of file GIDI_reaction.cc.

566 {
567
568 Vector vector( 0 );
569
570 if( !m_isPairProduction ) vector += m_outputChannel->multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, a_productID );
571
572 return( vector );
573}

Referenced by GIDI::ProtareSingle::multiGroupAverageMomentum(), and multiGroupDepositionMomentum().

◆ multiGroupCrossSection()

Vector GIDI::Reaction::multiGroupCrossSection ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_label = "" ) const

Returns the multi-group, cross section for the requested label the reaction.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_label[in] If not an empty string, this is used as the label for the form to return and the a_temperatureInfo labels are ignored.
Returns
The requested multi-group cross section as a GIDI::Vector.

Definition at line 349 of file GIDI_reaction.cc.

350 {
351
352 Vector vector( 0 );
353
354 Functions::Gridded1d const *form = dynamic_cast<Functions::Gridded1d const*>(
355 a_settings.form( a_smr, m_crossSection, a_temperatureInfo, "cross section", a_label ) );
356 if( form != nullptr ) vector = form->data( );
357
358 return( vector );
359}

Referenced by multiGroupAverageEnergy(), GIDI::ProtareSingle::multiGroupCrossSection(), multiGroupGain(), multiGroupMultiplicity(), and multiGroupProductMatrix().

◆ multiGroupDepositionEnergy()

Vector GIDI::Reaction::multiGroupDepositionEnergy ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
Transporting::Particles const & a_particles ) const

Returns the multi-group, total deposition energy for the requested label for this reaction. This is a cross section weighted deposition energy. The deposition energy is calculated by subtracting the average energy from each transportable particle from the available energy. The list of transportable particles is specified via the list of particle specified in the a_settings argument. This method does not include any photon deposition energy for this reaction that is in the GNDS orphanProducts node.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label and the products that are transported.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_particles[in] The list of particles to be transported.
Returns
The requested multi-group deposition energy as a GIDI::Vector.

Definition at line 509 of file GIDI_reaction.cc.

510 {
511
512 std::map<std::string, Transporting::Particle> const &products = a_particles.particles( );
513 Vector vector;
514
515 if( moniker( ) != GIDI_orphanProductChars ) {
516 if( ( a_settings.zeroDepositionIfAllProductsTracked( ) ) && areAllProductsTracked( a_particles ) && !m_isPairProduction ) return( vector );
517
518 vector = multiGroupAvailableEnergy( a_smr, a_settings, a_temperatureInfo );
519 }
520
521 Vector availableEnergy( vector );
522
523 for( std::map<std::string, Transporting::Particle>::const_iterator iter = products.begin( ); iter != products.end( ); ++iter ) {
524 vector -= multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, iter->first );
525 }
526
527 for( std::size_t index = 0; index < availableEnergy.size( ); ++index ) { // Check for values that should probably be 0.0.
528 if( std::fabs( vector[index] ) < std::fabs( 1e-14 * availableEnergy[index] ) ) vector[index] = 0.0;
529 }
530
531 return( vector );
532}
#define GIDI_orphanProductChars
Definition GIDI.hpp:182
bool areAllProductsTracked(Transporting::Particles const &a_particles) const
Vector multiGroupAvailableEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
Vector multiGroupAverageEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Component & availableEnergy()
Definition GIDI.hpp:4302
std::string const & moniker() const
Definition GUPI.hpp:102

Referenced by GIDI::ProtareSingle::multiGroupDepositionEnergy().

◆ multiGroupDepositionMomentum()

Vector GIDI::Reaction::multiGroupDepositionMomentum ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
Transporting::Particles const & a_particles ) const

Returns the multi-group, total deposition momentum for the requested label for this reaction. This is a cross section weighted deposition momentum. The deposition momentum is calculated by subtracting the average momentum from each transportable particle from the available momentum. The list of transportable particles is specified via the list of particle specified in the a_settings argument. This method does not include any photon deposition momentum for this reaction that is in the GNDS orphanProducts node.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_particles[in] The list of particles to be transported.
Returns
The requested multi-group deposition momentum as a GIDI::Vector.

Definition at line 589 of file GIDI_reaction.cc.

590 {
591
592 std::map<std::string, Transporting::Particle> const &products = a_particles.particles( );
593 Vector vector;
594
595 if( moniker( ) != GIDI_orphanProductChars ) {
596 vector = multiGroupAvailableMomentum( a_smr, a_settings, a_temperatureInfo );
597 }
598
599 for( std::map<std::string, Transporting::Particle>::const_iterator iter = products.begin( ); iter != products.end( ); ++iter ) {
600 vector -= multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, iter->first );
601 }
602
603 return( vector );
604}
Vector multiGroupAvailableMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
Vector multiGroupAverageMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const

◆ multiGroupFissionMatrix()

Matrix GIDI::Reaction::multiGroupFissionMatrix ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
Transporting::Particles const & a_particles,
std::size_t a_order ) const

Like Reaction::multiGroupProductMatrix, but only returns the fission neutron, transfer matrix.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label and if delayed neutrons should be included.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_particles[in] The list of particles to be transported.
a_order[in] Requested product matrix, Legendre order.
Returns
The requested multi-group neutron fission matrix as a GIDI::Matrix.

Definition at line 436 of file GIDI_reaction.cc.

437 {
438
439 Matrix matrix( 0, 0 );
440
441 if( hasFission( ) ) matrix += multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, PoPI::IDs::neutron, a_order );
442 return( matrix );
443}
Matrix multiGroupProductMatrix(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::string const &a_productID, std::size_t a_order) const
static std::string const neutron
Definition PoPI.hpp:164

Referenced by GIDI::ProtareSingle::multiGroupFissionMatrix().

◆ multiGroupGain()

Vector GIDI::Reaction::multiGroupGain ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_productID,
std::string const & a_projectileID ) const

Returns the multi-group, gain for the requested particle and label. This is a cross section weighted gain summed over all reactions. If a_productID and a_projectileID are the same, then the multi-group cross section is subtracted for the returned value to indicate that the a_projectileID as been absorted.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_productID[in] The particle PoPs' id for the whose gain is to be calculated.
a_projectileID[in] The particle PoPs' id for the projectile.
Returns
The requested multi-group gain as a GIDI::Vector.

Definition at line 620 of file GIDI_reaction.cc.

621 {
622
623 Vector vector( multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID ) );
624
625 if( PoPI::compareSpecialParticleIDs( a_productID, a_projectileID ) ) vector -= multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo );
626
627 return( vector );
628}
Vector multiGroupMultiplicity(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
bool compareSpecialParticleIDs(std::string const &a_id1, std::string const &a_id2)
Definition PoPI_misc.cc:133

Referenced by GIDI::ProtareSingle::multiGroupGain().

◆ multiGroupMultiplicity()

Vector GIDI::Reaction::multiGroupMultiplicity ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
std::string const & a_productID ) const

Returns the multi-group, total multiplicity for the requested label for the requested product. This is a cross section weighted multiplicity.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_productID[in] Particle id for the requested product.
Returns
The requested multi-group multiplicity as a GIDI::Vector.

Definition at line 226 of file GIDI_reaction.cc.

227 {
228
229 Vector vector( 0 );
230
231 if( m_decayPositronium ) {
232 if( a_productID == PoPI::IDs::photon ) vector += multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo ) * 2; }
233 else {
234 vector += m_outputChannel->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID );
235 }
236
237 return( vector );
238}

Referenced by GIDI::ProtareSingle::multiGroupFissionGammaMultiplicity(), GIDI::ProtareSingle::multiGroupFissionNeutronMultiplicity(), multiGroupGain(), and GIDI::ProtareSingle::multiGroupMultiplicity().

◆ multiGroupProductMatrix()

Matrix GIDI::Reaction::multiGroupProductMatrix ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
Transporting::Particles const & a_particles,
std::string const & a_productID,
std::size_t a_order ) const

Returns the multi-group, product matrix for the requested label for the requested product index for the requested Legendre order. If no data are found, an empty GIDI::Matrix is returned.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label and if delayed neutrons should be included.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_particles[in] The list of particles to be transported.
a_productID[in] Particle id for the requested product.
a_order[in] Requested product matrix, Legendre order.
Returns
The requested multi-group product matrix as a GIDI::Matrix.

Definition at line 395 of file GIDI_reaction.cc.

397 {
398
399 Matrix matrix( 0, 0 );
400
401 if( m_decayPositronium ) {
402 if( a_productID == PoPI::IDs::photon ) {
403 if( a_order == 0 ) {
404 Vector productionCrossSection = multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo ) * 2;
405 std::map<std::string, GIDI::Transporting::Particle> const &particles = a_particles.particles( );
406 std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = particles.find( PoPI::IDs::photon );
407 GIDI::Transporting::MultiGroup const &multiGroup = particle->second.fineMultiGroup( );
408 int multiGroupIndexFromEnergy = multiGroup.multiGroupIndexFromEnergy( PoPI_electronMass_MeV_c2, true );
409 Matrix matrix2( productionCrossSection.size( ), productionCrossSection.size( ) );
410
411 for( std::size_t i1 = 0; i1 < productionCrossSection.size( ); ++i1 ) {
412 matrix2.set( i1, static_cast<std::size_t>( multiGroupIndexFromEnergy ), productionCrossSection[i1] );
413 }
414 matrix += matrix2;
415 }
416 } }
417 else {
418 matrix += m_outputChannel->multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_productID, a_order );
419 }
420
421 return( matrix );
422}

Referenced by multiGroupFissionMatrix(), and GIDI::ProtareSingle::multiGroupProductMatrix().

◆ multiGroupQ()

Vector GIDI::Reaction::multiGroupQ ( LUPI::StatusMessageReporting & a_smr,
Transporting::MG const & a_settings,
Styles::TemperatureInfo const & a_temperatureInfo,
bool a_final ) const

Returns the multi-group, total Q for the requested label. This is a cross section weighted Q.

Parameters
a_smr[Out] If errors are not to be thrown, then the error is reported via this instance.
a_settings[in] Specifies the requested label.
a_temperatureInfo[in] Specifies the temperature and labels use to lookup the requested data.
a_final[in] If false, only the Q for the primary reactions are return, otherwise, the Q for the final reactions.
a_reactionsToExclude[in] A list of reaction indices that are to be ignored when calculating the Q.
Returns
The requested multi-group Q as a GIDI::Vector.

Definition at line 373 of file GIDI_reaction.cc.

374 {
375
376 if( m_decayPositronium ) return( Vector { 0 } ); // Special case, returns Q with all 0.0s.
377
378 return( m_outputChannel->multiGroupQ( a_smr, a_settings, a_temperatureInfo, a_final ) );
379}

Referenced by GIDI::ProtareSingle::multiGroupQ().

◆ nuclearPlusInterferencePresent()

bool GIDI::Reaction::nuclearPlusInterferencePresent ( ) const
inline

Returns the value of m_nuclearPlusInterferencePresent member.

Definition at line 4294 of file GIDI.hpp.

◆ onlyRutherfordScatteringPresent()

bool GIDI::Reaction::onlyRutherfordScatteringPresent ( ) const
inline

Returns the value of m_onlyRutherfordScatteringPresent member.

Definition at line 4292 of file GIDI.hpp.

◆ outputChannel()

OutputChannel * GIDI::Reaction::outputChannel ( ) const
inline

Returns a reference to the m_outputChannel member.

Definition at line 4307 of file GIDI.hpp.

Referenced by GIDI::ProtareSingle::isDelayedFissionNeutronComplete().

◆ productIDs()

void GIDI::Reaction::productIDs ( std::set< std::string > & a_ids,
Transporting::Particles const & a_particles,
bool a_transportablesOnly ) const

Fills in a std::set with a unique list of all product indices produced by this reaction. If a_transportablesOnly is true, only transportable product indices are return.

Parameters
a_ids[out] The unique list of product indices.
a_particles[in] The list of particles to be transported.
a_transportablesOnly[in] If true, only transportable product indices are added in the list.

Definition at line 180 of file GIDI_reaction.cc.

180 {
181
182 if( m_decayPositronium ) {
183 std::string electronAnti( PoPI::IDs::electron + PoPI::IDs::anti );
184 std::set<std::string> ids;
185 m_outputChannel->productIDs( ids, a_particles, a_transportablesOnly );
186 for( auto iterId = ids.begin( ); iterId != ids.end( ); ++iterId ) {
187 if( ( *iterId == PoPI::IDs::electron ) || ( *iterId == electronAnti ) ) continue;
188 a_ids.insert( *iterId );
189 }
190 a_ids.insert( PoPI::IDs::photon ); }
191 else {
192 m_outputChannel->productIDs( a_ids, a_particles, a_transportablesOnly );
193 }
194}
static std::string const anti
Definition PoPI.hpp:173
static std::string const electron
Definition PoPI.hpp:163

Referenced by GIDI::ProtareSingle::productIDs().

◆ productMultiplicity()

int GIDI::Reaction::productMultiplicity ( std::string const & a_productID) const
inline

Returns the product multiplicity (e.g., 0, 1, 2, ...) or -1 if energy dependent or not an integer.

Definition at line 4318 of file GIDI.hpp.

◆ reactionIndex()

std::size_t GIDI::Reaction::reactionIndex ( ) const
inline

Returns the value of the m_reactionIndex member.

Definition at line 4282 of file GIDI.hpp.

◆ recalculateMultiGroupData()

void GIDI::Reaction::recalculateMultiGroupData ( ProtareSingle const * a_protare,
Styles::TemperatureInfo const & a_temperatureInfo )

This function recalculates the multi-group data for data labelled with a_temperatureInfo.heatedMultiGroup().

Parameters
a_temperatureInfo[in] The temperature for which multi-group data are to be recalculated.

Definition at line 831 of file GIDI_reaction.cc.

831 {
832
833 std::string heatedMultiGroupLabel = a_temperatureInfo.heatedMultiGroup( );
834
835 GUPI::Ancestry const *ancestorRoot = root( );
836 if( ancestorRoot->moniker( ) != GIDI_topLevelChars ) throw Exception( "Reaction::recalculateMultiGroupData: could not find parent protare." );
837
838 Styles::Suite const &styles = a_protare->styles( );
839 if( styles.find( heatedMultiGroupLabel ) == styles.end( ) ) return;
840
841 Styles::HeatedMultiGroup const &heatedMultiGroupStyle = *styles.get<Styles::HeatedMultiGroup const>( heatedMultiGroupLabel );
842
843 std::vector<double> groupBoundaries = heatedMultiGroupStyle.groupBoundaries( a_protare->projectile( ).ID( ) );
844 Transporting::MultiGroup multiGroup = Transporting::MultiGroup( "recal", groupBoundaries );
845
846 Transporting::Flux flux( "recal", a_temperatureInfo.temperature( ).value( ) );
847 double energies[2] = { groupBoundaries[0], groupBoundaries.back( ) };
848 double fluxValues[2] = { 1.0, 1.0 };
849 Transporting::Flux_order flux_order( 0, 2, energies, fluxValues );
850 flux.addFluxOrder( flux_order );
851
852 MultiGroupCalulationInformation multiGroupCalulationInformation( multiGroup, flux );
853 calculateMultiGroupData( a_protare, a_temperatureInfo, heatedMultiGroupLabel, multiGroupCalulationInformation );
854}
#define GIDI_topLevelChars
Definition GIDI.hpp:169
void calculateMultiGroupData(ProtareSingle const *a_protare, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_heatedMultiGroupLabel, MultiGroupCalulationInformation const &a_multiGroupCalulationInformation)

Referenced by modifyCrossSection().

◆ RutherfordScatteringPresent()

bool GIDI::Reaction::RutherfordScatteringPresent ( ) const
inline

Returns the value of m_RutherfordScatteringPresent member.

Definition at line 4290 of file GIDI.hpp.

Referenced by GIDI::ProtareSingle::checkIf_nuclearPlusCoulombInterferenceWanted().

◆ setActive()

void GIDI::Reaction::setActive ( bool a_active)
inline

Sets m_active to a_active.

Definition at line 4281 of file GIDI.hpp.

◆ setOutputChannel()

void GIDI::Reaction::setOutputChannel ( OutputChannel * a_outputChannel)

Sets this reaction's output channel to a_outputChannel.

Parameters
a_outputChannel[in] The output channel to make this reaction output channel.

Definition at line 263 of file GIDI_reaction.cc.

263 {
264
265 m_outputChannel = a_outputChannel;
266 m_outputChannel->setAncestor( this );
267}

◆ threshold()

double GIDI::Reaction::threshold ( ) const
inline

Returns the value of the m_QThreshold member.

Definition at line 4323 of file GIDI.hpp.

◆ toXMLList()

void GIDI::Reaction::toXMLList ( GUPI::WriteInfo & a_writeInfo,
std::string const & a_indent = "" ) const
virtual

Fills the argument a_writeInfo with the XML lines that represent this. Recursively enters each sub-node.

Parameters
a_writeInfo[in/out] Instance containing incremental indentation and other information and stores the appended lines.
a_indent[in] The amount to indent this node.

Reimplemented from GUPI::Ancestry.

Definition at line 894 of file GIDI_reaction.cc.

894 {
895
896 std::string indent2 = a_writeInfo.incrementalIndent( a_indent );
897 std::string attributes;
898
899 attributes += a_writeInfo.addAttribute( GIDI_labelChars, label( ) );
900 attributes += a_writeInfo.addAttribute( GIDI_ENDF_MT_Chars, intToString( ENDF_MT( ) ) );
901 if( m_fissionGenre != "" ) attributes += a_writeInfo.addAttribute( GIDI_fissionGenreChars, m_fissionGenre );
902 a_writeInfo.addNodeStarter( a_indent, moniker( ), attributes );
903
904 m_doubleDifferentialCrossSection.toXMLList( a_writeInfo, indent2 );
905 m_crossSection.toXMLList( a_writeInfo, indent2 );
906 if( m_outputChannel != nullptr ) m_outputChannel->toXMLList( a_writeInfo, indent2 );
907 m_availableEnergy.toXMLList( a_writeInfo, indent2 );
908 m_availableMomentum.toXMLList( a_writeInfo, indent2 );
909
910 a_writeInfo.addNodeEnder( moniker( ) );
911}
int ENDF_MT() const
Definition GIDI.hpp:4284
void addNodeEnder(std::string const &a_moniker)
Definition GUPI.hpp:59
std::string incrementalIndent(std::string const &indent)
Definition GUPI.hpp:52
void addNodeStarter(std::string const &indent, std::string const &a_moniker, std::string const &a_attributes="")
Definition GUPI.hpp:55
std::string addAttribute(std::string const &a_name, std::string const &a_value) const
Definition GUPI.hpp:60
std::string intToString(int a_value)
Definition GIDI_misc.cc:415

◆ twoBodyThreshold()

double GIDI::Reaction::twoBodyThreshold ( ) const
inline

Returns the value of the m_twoBodyThreshold member.

Definition at line 4325 of file GIDI.hpp.

◆ xlinkItemKey()

std::string GIDI::Reaction::xlinkItemKey ( ) const
inlinevirtual

Returns the result of calling "GUPI::Ancestry::buildXLinkItemKey( GIDI_labelChars, label() )".

Reimplemented from GIDI::Form.

Definition at line 4314 of file GIDI.hpp.

◆ ProtareSingle

friend class ProtareSingle
friend

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