Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
GIDI_protare.cc
Go to the documentation of this file.
1/*
2# <<BEGIN-copyright>>
3# Copyright 2019, Lawrence Livermore National Security, LLC.
4# This file is part of the gidiplus package (https://github.com/LLNL/gidiplus).
5# gidiplus is licensed under the MIT license (see https://opensource.org/licenses/MIT).
6# SPDX-License-Identifier: MIT
7# <<END-copyright>>
8*/
9
10#include <stdlib.h>
11#include <algorithm>
12#include <cmath>
13
14#include "GIDI.hpp"
15#include <HAPI.hpp>
16
17namespace GIDI {
18
19static bool sortTemperatures( Styles::TemperatureInfo const &lhs, Styles::TemperatureInfo const &rhs );
20
21/*! \class Protare
22 * Base class for the protare sub-classes.
23 */
24
25/* *********************************************************************************************************//**
26 * Base Protare constructor.
27 ***********************************************************************************************************/
28
30 GUPI::Ancestry( "" ),
31 m_projectile( "", "", -1.0 ),
32 m_target( "", "", -1.0 ),
33 m_GNDS_target( "", "", -1.0 ) {
34
35}
36
37/* *********************************************************************************************************//**
38 ***********************************************************************************************************/
39
43
44/* *********************************************************************************************************//**
45 * Called by the constructs. This method does most of the parsing.
46 *
47 * @param a_node [in] The protare (i.e., reactionSuite) node to be parsed and used to construct a Protare.
48 * @param a_setupInfo [in] Information create my the Protare constructor to help in parsing.
49 * @param a_pops [in] A PoPs Database instance used to get particle indices and possibly other particle information.
50 * @param a_internalPoPs [in] The internal PoPI::Database instance used to get particle indices and possibly other particle information.
51 * @param a_targetRequiredInGlobalPoPs [in] If *true*, the target is required to be in **a_pops**.
52 * @param a_requiredInPoPs [in] If *true*, particle is required to be in **a_pops**.
53 ***********************************************************************************************************/
54
55void Protare::initialize( HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pops, PoPI::Database const &a_internalPoPs,
56 bool a_targetRequiredInGlobalPoPs, bool a_requiredInPoPs ) {
57
58 setMoniker( a_node.name( ) );
59
60 std::string projectileID = a_node.attribute_as_string( GIDI_projectileChars );
61 m_projectile = ParticleInfo( projectileID, a_pops, a_internalPoPs, a_requiredInPoPs );
62
63 std::string targetID = a_node.attribute_as_string( GIDI_targetChars );
64 m_GNDS_target = ParticleInfo( targetID, a_pops, a_internalPoPs, a_targetRequiredInGlobalPoPs && a_requiredInPoPs );
65
66 auto iter = a_setupInfo.m_particleSubstitution->find( m_GNDS_target.ID( ) );
67 if( iter != a_setupInfo.m_particleSubstitution->end( ) ) {
68 m_target = iter->second; }
69 else {
70 m_target = m_GNDS_target;
71 }
72}
73
74/* *********************************************************************************************************//**
75 * If the protare is a ProtareTNSL then summing over all reactions will include the standard protare's elastic cross section
76 * in the domain of the TNSL data. The standard elastic cross section should not be added in this domain.
77 * If needed, this function corrects the cross section for this over counting of the elastic cross section.
78 * This method does nothing unless overwritten by the ProtareTNSL class.
79 *
80 * @param a_label [in] The label of the elastic cross section data to use if over counting needs to be corrected.
81 * @param a_crossSectionSum [in] The cross section to correct.
82 ***********************************************************************************************************/
83
84void Protare::TNSL_crossSectionSumCorrection( LUPI_maybeUnused std::string const &a_label, LUPI_maybeUnused Functions::XYs1d &a_crossSectionSum ) {
85
86}
87
88/* *********************************************************************************************************//**
89 * If the protare is a ProtareTNSL then summing over all reactions will include the standard protare's elastic cross section
90 * in the domain of the TNSL data. The standard elastic cross section should not be added in this domain.
91 * If needed, this function corrects the cross section for this over counting of the elastic cross section.
92 * This method does nothing unless overwritten by the ProtareTNSL class.
93 *
94 * @param a_label [in] The label of the elastic cross section data to use if over counting needs to be corrected.
95 * @param a_crossSectionSum [in] The cross section to correct.
96 ***********************************************************************************************************/
97
98void Protare::TNSL_crossSectionSumCorrection( LUPI_maybeUnused std::string const &a_label, LUPI_maybeUnused Functions::Ys1d &a_crossSectionSum ) {
99
100}
101
102/* *********************************************************************************************************//**
103 * If the protare is a ProtareTNSL then summing over all reactions will include the standard protare's elastic cross section
104 * in the domain of the TNSL data. The standard elastic cross section should not be added in this domain.
105 * If needed, this function corrects the cross section for this over counting of the elastic cross section.
106 * This method does nothing as the multi-group cross section data for the standard protare's elastic cross section are
107 * zeroed when the data are read in. However, this method is added so codes do not have to check the type of data they are accessing.
108 *
109 * @param a_label [in] The label of the elastic cross section data to use if over counting needs to be corrected.
110 * @param a_crossSectionSum [in] The cross section to correct.
111 ***********************************************************************************************************/
112
113void Protare::TNSL_crossSectionSumCorrection( LUPI_maybeUnused std::string const &a_label, LUPI_maybeUnused Vector &a_crossSectionSum ) {
114
115}
116
117/* *********************************************************************************************************//**
118 * Returns a list of all reaction indices whose ENDL C value is in the set *a_CValues*.
119 *
120 * @param a_CValues [in] A list of ENDL C values.
121 * @param a_checkActiveState [in] If true, all reactions whose active state is false are not included in the returned set even if their CValue match on in the list.
122 *
123 * @return The list of reaction indices.
124 ***********************************************************************************************************/
125
126ExcludeReactionsSet Protare::reactionIndicesMatchingENDLCValues( std::set<int> const &a_CValues, bool a_checkActiveState ) {
127
128 ExcludeReactionsSet indices;
129
130 for( std::size_t i1 = 0; i1 < numberOfReactions( ); ++i1 ) {
131 Reaction *reaction1 = reaction( i1 );
132
133 if( a_checkActiveState && !reaction1->active( ) ) continue;
134 if( a_CValues.find( reaction1->ENDL_C( ) ) != a_CValues.end( ) ) indices.insert( i1 );
135 }
136
137 return( indices );
138}
139
140/*! \class ProtareSingle
141 * Class to store a GNDS <**reactionSuite**> node.
142 */
143
144/* *********************************************************************************************************//**
145 * Parses a GNDS file to construct the Protare instance. Calls the initialize method which does most of the work.
146 *
147 * @param a_pops [in] A PoPs Database instance used to get particle indices and possibly other particle information.
148 * @param a_projectileID [in] The PoPs id of the projectile.
149 * @param a_targetID [in] The PoPs id of the target.
150 * @param a_evaluation [in] The evaluation string.
151 * @param a_interaction [in] The interaction flag for the protare.
152 * @param a_formatVersion [in] The GNDS format to use.
153 ***********************************************************************************************************/
154
155ProtareSingle::ProtareSingle( PoPI::Database const &a_pops, std::string const &a_projectileID, std::string const &a_targetID, std::string const &a_evaluation,
156 std::string const &a_interaction, std::string const &a_formatVersion ) :
157 m_doc( nullptr ),
158 m_dataManager( nullptr ),
159 m_numberOfLazyParsingHelperForms( 0 ),
160 m_numberOfLazyParsingHelperFormsReplaced( 0 ),
161 m_formatVersion( a_formatVersion ),
162 m_evaluation( a_evaluation ),
163 m_interaction( a_interaction ),
164 m_projectileFrame( Frame::lab ),
165 m_decayPositronium( false ),
166 m_thresholdFactor( 0.0 ),
167 m_nuclearPlusCoulombInterferenceOnlyReaction( nullptr ),
168 m_pointwiseAverageProductEnergy( GIDI_averageEnergyChars, GIDI_labelChars ) {
169
171 initialize( );
172
173 setProjectile( ParticleInfo( a_projectileID, a_pops, a_pops, true ) );
174 setTarget( ParticleInfo( a_targetID, a_pops, a_pops, true ) );
175}
176
177/* *********************************************************************************************************//**
178 * Parses a GNDS file to construct the Protare instance. Calls the initialize method which does most of the work.
179 *
180 * @param a_construction [in] Used to pass user options to the constructor.
181 * @param a_fileName [in] File containing a protare (i.e., reactionSuite) node that is parsed and used to construct the Protare.
182 * @param a_fileType [in] File type of a_fileType. Currently, only GIDI::FileType::XML and GIDI::FileType::HDF are supported.
183 * @param a_pops [in] A PoPs Database instance used to get particle indices and possibly other particle information.
184 * @param a_particleSubstitution [in] Map of particles to substitute with another particles.
185 * @param a_libraries [in] The list of libraries that were searched to find *this*.
186 * @param a_interaction [in] The interaction flag for the protare.
187 * @param a_targetRequiredInGlobalPoPs [in] If *true*, the target is required to be in **a_pops**.
188 * @param a_requiredInPoPs [in] If *true*, no particle is required to be in **a_pops**.
189 ***********************************************************************************************************/
190
191ProtareSingle::ProtareSingle( Construction::Settings const &a_construction, std::string const &a_fileName, FileType a_fileType,
192 PoPI::Database const &a_pops, ParticleSubstitution const &a_particleSubstitution, std::vector<std::string> const &a_libraries,
193 std::string const &a_interaction, bool a_targetRequiredInGlobalPoPs, bool a_requiredInPoPs ) :
194 Protare( ),
195 m_doc( nullptr ),
196 m_dataManager( nullptr ),
197 m_numberOfLazyParsingHelperForms( 0 ),
198 m_numberOfLazyParsingHelperFormsReplaced( 0 ),
199 m_libraries( a_libraries ),
200 m_interaction( a_interaction ),
201 m_fileName( a_fileName ),
202 m_realFileName( LUPI::FileInfo::realPath( a_fileName ) ),
203 m_decayPositronium( a_construction.decayPositronium( ) ),
204 m_pointwiseAverageProductEnergy( GIDI_averageEnergyChars, GIDI_labelChars ) {
205
206#ifdef HAPI_USE_PUGIXML
207 if( a_fileType == GIDI::FileType::XML ) {
208 m_doc = new HAPI::PugiXMLFile( a_fileName.c_str( ), "ProtareSingle::ProtareSingle" );
209 }
210#endif
211#ifdef HAPI_USE_HDF5
212 if( a_fileType == GIDI::FileType::HDF ) {
213 m_doc = new HAPI::HDFFile( a_fileName.c_str( ) );
214 }
215#endif
216 if( m_doc == nullptr ) {
217 throw std::runtime_error( "Only XML/HDF file types supported." );
218 }
219
220#ifdef GIDIP_TEST_PARSING
221 if( a_construction.parseMode( ) != Construction::ParseMode::noParsing ) {
222#endif
223 HAPI::Node protare = m_doc->first_child( );
224
225 SetupInfo setupInfo( this );
226 ParticleSubstitution particleSubstitution( a_particleSubstitution );
227 setupInfo.m_particleSubstitution = &particleSubstitution;
228
229 initialize( a_construction, protare, setupInfo, a_pops, a_targetRequiredInGlobalPoPs, a_requiredInPoPs );
230#ifdef GIDIP_TEST_PARSING
231 }
232#endif
233}
234
235/* *********************************************************************************************************//**
236 * Parses a GNDS HAPI::Node instance to construct the Protare instance. Calls the ProtareSingle::initialize method which does most of the work.
237 *
238 * @param a_construction [in] Used to pass user options to the constructor.
239 * @param a_node [in] The **HAPI::Node** to be parsed and used to construct the Protare.
240 * @param a_pops [in] A PoPI::Database instance used to get particle indices and possibly other particle information.
241 * @param a_particleSubstitution [in] Map of particles to substitute with another particles.
242 * @param a_libraries [in] The list of libraries that were searched to find *this*.
243 * @param a_targetRequiredInGlobalPoPs [in] If *true*, the target is required to be in **a_pops**.
244 * @param a_interaction [in] The interaction between the projectile and target.
245 * @param a_requiredInPoPs [in] If *true*, no particle is required to be in **a_pops**.
246 ***********************************************************************************************************/
247
248ProtareSingle::ProtareSingle( Construction::Settings const &a_construction, HAPI::Node const &a_node, PoPI::Database const &a_pops,
249 ParticleSubstitution const &a_particleSubstitution, std::vector<std::string> const &a_libraries,
250 LUPI_maybeUnused std::string const &a_interaction, bool a_targetRequiredInGlobalPoPs, bool a_requiredInPoPs ) :
251 Protare( ),
252 m_doc( nullptr ),
253 m_dataManager( nullptr ),
254 m_numberOfLazyParsingHelperForms( 0 ),
255 m_numberOfLazyParsingHelperFormsReplaced( 0 ),
256 m_libraries( a_libraries ),
257 m_pointwiseAverageProductEnergy( GIDI_averageEnergyChars, GIDI_labelChars ) {
258
259 SetupInfo setupInfo( this );
260 ParticleSubstitution particleSubstitution( a_particleSubstitution );
261 setupInfo.m_particleSubstitution = &particleSubstitution;
262
263 initialize( a_construction, a_node, setupInfo, a_pops, a_targetRequiredInGlobalPoPs, a_requiredInPoPs );
264}
265
266/* *********************************************************************************************************//**
267 * Base initializer to be called by all constructors (directly or indirectly).
268 ***********************************************************************************************************/
269
270void ProtareSingle::initialize( ) {
271
272 m_externalFiles.setMoniker( GIDI_externalFilesChars );
273 m_externalFiles.setAncestor( this );
274
275 m_styles.setMoniker( GIDI_stylesChars );
276 m_styles.setAncestor( this );
277
278 m_documentations.setAncestor( this );
279 m_documentations.setMoniker( GIDI_documentations_1_10_Chars );
280
281 m_reactions.setAncestor( this );
282 m_reactions.setMoniker( GIDI_reactionsChars );
283
284 m_orphanProducts.setAncestor( this );
285 m_orphanProducts.setMoniker( GIDI_orphanProductsChars );
286
287 m_incompleteReactions.setAncestor( this );
288 m_incompleteReactions.setMoniker( GIDI_reactionsChars );
289
290 m_sums.setAncestor( this );
291
292 m_fissionComponents.setAncestor( this );
293 m_fissionComponents.setMoniker( GIDI_fissionComponentsChars );
294
295 m_RutherfordScatteringPresent = false;
296 m_onlyRutherfordScatteringPresent = false;
297 m_nuclearPlusCoulombInterferenceOnlyReaction = nullptr;
298 m_multiGroupSummedReaction = nullptr;
299 m_multiGroupSummedDelayedNeutrons = nullptr;
300
301 m_ACE_URR_probabilityTables.setAncestor( this );
302 m_ACE_URR_probabilityTables.setMoniker( GIDI_ACE_URR_probabilityTablesChars );
303
304 m_photoAtomicIncoherentDoppler.setAncestor( this );
305 m_photoAtomicIncoherentDoppler.setMoniker( GIDI_LLNL_photoAtomicIncoherentDoppler_Chars );
306
307 m_pointwiseAverageProductEnergy.setAncestor( this );
308 m_GRIN_continuumGammas = nullptr;
309}
310
311/* *********************************************************************************************************//**
312 * Called by the constructs. This method does most of the parsing.
313 *
314 * @param a_construction [in] Used to pass user options to the constructor.
315 * @param a_node [in] The protare (i.e., reactionSuite) node to be parsed and used to construct a Protare.
316 * @param a_setupInfo [in] Information create my the Protare constructor to help in parsing.
317 * @param a_pops [in] A PoPs Database instance used to get particle indices and possibly other particle information.
318 * @param a_targetRequiredInGlobalPoPs [in] If *true*, the target is required to be in **a_pops**.
319 * @param a_requiredInPoPs [in] If *true*, no particle is required to be in **a_pops**.
320 ***********************************************************************************************************/
321
322void ProtareSingle::initialize( Construction::Settings const &a_construction, HAPI::Node const &a_node, SetupInfo &a_setupInfo,
323 PoPI::Database const &a_pops, bool a_targetRequiredInGlobalPoPs, bool a_requiredInPoPs ) {
324
325 if( a_node.name( ) != GIDI_topLevelChars ) throw Exception( "Node '" + a_node.name( ) + "' is not a 'reactionSuite' node." );
326
327 m_externalFiles.parse( a_construction, a_node.child( GIDI_externalFilesChars ), a_setupInfo, a_pops, m_internalPoPs, parseExternalFilesSuite, nullptr );
328 std::string parentDir = m_fileName.substr( 0, m_fileName.find_last_of( "/" ) );
329 m_externalFiles.registerBinaryFiles( parentDir, a_setupInfo );
330
331 HAPI::Node internalPoPs = a_node.child( GIDI_PoPsChars );
332 m_internalPoPs.addDatabase( internalPoPs, true );
333 std::vector<PoPI::Alias *> const aliases = m_internalPoPs.aliases( );
334 for( auto alias = aliases.begin( ); alias != aliases.end( ); ++alias ) {
335 a_setupInfo.m_particleSubstitution->insert( { (*alias)->pid( ), ParticleInfo( (*alias)->ID( ), a_pops, a_pops, true ) } );
336 }
337
338 m_interaction = a_node.attribute_as_string( GIDI_interactionChars );
339 if( m_interaction == GIDI_MapInteractionTNSLChars ) a_targetRequiredInGlobalPoPs = false;
340
341 Protare::initialize( a_node, a_setupInfo, a_pops, m_internalPoPs, a_targetRequiredInGlobalPoPs, a_requiredInPoPs );
342 initialize( );
343
344 m_formatVersion.setFormat( a_node.attribute_as_string( GIDI_formatChars ) );
345 if( !m_formatVersion.supported( ) ) throw Exception( "unsupported GNDS format version" );
346 a_setupInfo.m_formatVersion = m_formatVersion;
347
348 if( m_formatVersion.major( ) > 1 ) {
349 m_interaction = a_node.attribute_as_string( GIDI_interactionChars ); }
350 else {
351 HAPI::Node const firstReaction = a_node.child( GIDI_reactionsChars ).first_child( );
352 if( firstReaction.attribute_as_int( GIDI_ENDF_MT_Chars ) == 502 ) m_interaction = GIDI_MapInteractionAtomicChars;
353 }
354 m_isPhotoAtomic = ( m_interaction == GIDI_MapInteractionAtomicChars ) && ( PoPI::IDs::photon == projectile( ).ID( ) );
355
356 PoPI::Database const *GRIN_pops = nullptr;
357 HAPI::Node const &applicationData = a_node.child( GIDI_applicationDataChars );
358 std::vector<std::string> extraGammaBranchStates;
359 for( HAPI::Node child1 = applicationData.first_child( ); !child1.empty( ); child1.to_next_sibling( ) ) {
360 std::string nodeName( child1.name( ) );
361 if( nodeName == GIDI_institutionChars ) {
362 std::string label = child1.attribute_as_string( GIDI_labelChars );
363 if( label == GIDI_LLNL_GRIN_continuumGammas ) {
364 HAPI::Node child2 = child1.child( GIDI_GRIN_continuumGammasChars );
365 if( a_construction.GRIN_continuumGammas( ) ) {
366 m_GRIN_continuumGammas = new GRIN::GRIN_continuumGammas( a_construction, child2, a_setupInfo, a_pops, m_internalPoPs, *this, &m_styles );
367 GRIN_pops = &(m_GRIN_continuumGammas->pops( ));
368 extraGammaBranchStates.push_back( m_GRIN_continuumGammas->captureResidualId( ) );
369 }
370 }
371 }
372 }
373 m_internalPoPs.calculateNuclideGammaBranchStateInfos( m_nuclideGammaBranchStateInfos, GRIN_pops, extraGammaBranchStates );
374
375 m_isTNSL_ProtareSingle = false;
376 if( m_interaction == GIDI_TNSLChars ) m_interaction = GIDI_MapInteractionTNSLChars;
377 if( m_interaction == GIDI_MapInteractionTNSLChars ) {
378 m_isTNSL_ProtareSingle = true; }
379 else { // For some legacy GNDS 1.10 files.
381 m_isTNSL_ProtareSingle = name.find( "thermalNeutronScatteringLaw" ) != std::string::npos;
382 if( m_isTNSL_ProtareSingle ) m_interaction = GIDI_MapInteractionTNSLChars;
383 }
384
385 m_thresholdFactor = 1.0;
386 if( a_pops.exists( target( ).pid( ) ) ) {
387 m_thresholdFactor = 1.0 + projectile( ).mass( "amu" ) / target( ).mass( "amu" ); // BRB FIXME, I think only this statement needs to be in this if section.
388 }
389
390 m_evaluation = a_node.attribute_as_string( GIDI_evaluationChars );
391
392 m_projectileFrame = parseFrame( a_node, a_setupInfo, GIDI_projectileFrameChars );
393
394 m_styles.parse( a_construction, a_node.child( GIDI_stylesChars ), a_setupInfo, a_pops, m_internalPoPs, parseStylesSuite, nullptr );
395 m_styles.updateChainEnds( );
396
397 Styles::Evaluated *evaluated = m_styles.get<Styles::Evaluated>( 0 );
398
399 m_projectileEnergyMin = evaluated->projectileEnergyDomain( ).minimum( );
400 m_projectileEnergyMax = evaluated->projectileEnergyDomain( ).maximum( );
401
402 m_reactions.parse( a_construction, a_node.child( GIDI_reactionsChars ), a_setupInfo, a_pops, m_internalPoPs, parseReaction, &m_styles );
403 m_orphanProducts.parse( a_construction, a_node.child( GIDI_orphanProductsChars ), a_setupInfo, a_pops, m_internalPoPs, parseOrphanProduct, &m_styles );
404 m_incompleteReactions.parse( a_construction, a_node.child( GIDI_incompleteReactionsChars ), a_setupInfo, a_pops, m_internalPoPs, parseReaction, &m_styles );
405
406 m_sums.parse( a_construction, a_node.child( GIDI_sumsChars ), a_setupInfo, a_pops, m_internalPoPs );
407 m_fissionComponents.parse( a_construction, a_node.child( GIDI_fissionComponentsChars ), a_setupInfo, a_pops, m_internalPoPs, parseFissionComponent, &m_styles );
408
409 m_RutherfordScatteringPresent = false;
410 m_onlyRutherfordScatteringPresent = false;
411 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
412 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
413 if( reaction1->RutherfordScatteringPresent( ) ) {
414 m_RutherfordScatteringPresent = true;
415 m_onlyRutherfordScatteringPresent = reaction1->onlyRutherfordScatteringPresent( );
416 break;
417 }
418 }
419
420 for( HAPI::Node child1 = applicationData.first_child( ); !child1.empty( ); child1.to_next_sibling( ) ) {
421 std::string nodeName( child1.name( ) );
422
423 if( nodeName == GIDI_institutionChars ) {
424 std::string label = child1.attribute_as_string( GIDI_labelChars );
425 if( label == GIDI_LLNL_Chars ) {
426 for( HAPI::Node child2 = child1.first_child( ); !child2.empty( ); child2.to_next_sibling( ) ) {
428 HAPI::Node const reactionNode = child2.child( GIDI_reactionChars );
429 a_setupInfo.m_isENDL_C_9 = true;
430 m_nuclearPlusCoulombInterferenceOnlyReaction = new Reaction( a_construction, reactionNode, a_setupInfo, a_pops, m_internalPoPs, *this, &m_styles );
431 a_setupInfo.m_isENDL_C_9 = false;
432 bool dropC_9 = false;
433 auto &crossSectionSuite = m_nuclearPlusCoulombInterferenceOnlyReaction->crossSection( );
434 for( auto crossSectionIter = crossSectionSuite.begin( ); crossSectionIter != crossSectionSuite.end( ); ++crossSectionIter ) {
435 auto iter = crossSectionSuite.checkLazyParsingHelperFormIterator( crossSectionIter );
436 if( (*iter)->type( ) == FormType::XYs1d ) {
437 Functions::XYs1d *xys1d = static_cast<Functions::XYs1d *>( *iter );
438 auto ys = xys1d->ys( );
439 for( auto yIter = ys.begin( ); yIter != ys.end( ); ++yIter ) {
440 if( *yIter < 0.0 ) {
441 dropC_9 = true;
442 break;
443 }
444 }
445 break;
446 }
447 }
448 if( dropC_9 ) {
449 delete m_nuclearPlusCoulombInterferenceOnlyReaction;
450 m_nuclearPlusCoulombInterferenceOnlyReaction = nullptr;
451 }
452 }
453 } }
454 else if( label == GIDI_LLNL_multiGroupReactions_Chars ) {
455 HAPI::Node child2 = child1.child( GIDI_reactionChars );
456 m_multiGroupSummedReaction = new Reaction( a_construction, child2, a_setupInfo, a_pops, m_internalPoPs, *this, &m_styles ); }
458 HAPI::Node child2 = child1.child( GIDI_outputChannelChars );
459 m_multiGroupSummedDelayedNeutrons = new OutputChannel( a_construction, child2, a_setupInfo, a_pops, m_internalPoPs, &m_styles, true, false ); }
460 else if( label == GIDI_LLNL_URR_probability_tables_Chars ) {
461 m_ACE_URR_probabilityTables.parse( a_construction, child1.child( GIDI_ACE_URR_probabilityTablesChars ), a_setupInfo, a_pops,
462 m_internalPoPs, parseACE_URR_probabilityTables, &m_styles ); }
464 HAPI::Node child2 = child1.child( GIDI_reactionsChars );
465 m_photoAtomicIncoherentDoppler.parse( a_construction, child2, a_setupInfo, a_pops, m_internalPoPs, parseReaction, &m_styles );
466 if( a_construction.usePhotoAtomicIncoherentDoppler( ) && m_photoAtomicIncoherentDoppler.size( ) > 0 ) {
467 // Add the consistent doppler broadened incoherent reactions to the list.
468 while( m_photoAtomicIncoherentDoppler.size( ) > 0 ) {
469 Reaction *photoAtomicIncoherentDopplerReaction = m_photoAtomicIncoherentDoppler.pop<Reaction>( 0 );
470 m_reactions.add( photoAtomicIncoherentDopplerReaction );
471 }
472 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) { // Set standard incoherent reaction to inactive
473 Reaction *reaction1 = m_reactions.get<Reaction>( i1 );
474 if( reaction1->ENDF_MT( ) == 504 ) {
475 reaction1->setActive(false);
476 }
477 }
478 }
479 }
481 m_pointwiseAverageProductEnergy.parse( a_construction, child1.child( GIDI_averageEnergyChars ), a_setupInfo, a_pops,
482 m_internalPoPs, parseAverageEnergySuite, &m_styles ); }
483 else if( label == GIDI_LLNL_GRIN_continuumGammas ) { // Already parsed above.
484 continue;
485 } }
486 else {
487 std::cout << "parseStylesSuite: Ignoring unsupported style = '" << nodeName << "'." << std::endl;
488 }
489 }
490}
491
492/* *********************************************************************************************************//**
493 ***********************************************************************************************************/
494
496
497 delete m_doc;
498 delete m_dataManager;
499 delete m_nuclearPlusCoulombInterferenceOnlyReaction;
500 delete m_multiGroupSummedReaction;
501 delete m_multiGroupSummedDelayedNeutrons;
502}
503
504/* *********************************************************************************************************//**
505 * Checks various things to determine if it is okay to use summed multi-group data or not.
506 *
507 * @param a_settings [in] Specifies user options.
508 *
509 * @return Returns **true** if summed data are to be returned and **false** otherwise.
510 ***********************************************************************************************************/
511
512bool ProtareSingle::useMultiGroupSummedData( Transporting::MG const &a_settings, ExcludeReactionsSet const &a_reactionsToExclude ) const {
513
514 if( ( numberOfInactiveReactions( ) > 0 ) || ( a_reactionsToExclude.size( ) > 0 ) ) return( false );
515 if( m_multiGroupSummedReaction == nullptr ) return( false );
516 if( a_settings.nuclearPlusCoulombInterferenceOnly( ) && m_RutherfordScatteringPresent ) return( false );
517 if( m_decayPositronium ) return( false ); // Need to handle decay of positronium by reaction.
518
519 return( a_settings.useMultiGroupSummedData( ) );
520}
521
522/* *********************************************************************************************************//**
523 * Returns **true** if delayed neutrons are to be included and **false** otherwise.
524 *
525 * @param a_settings [in] Specifies user options.
526 *
527 * @return Returns **true** if delayed neutrons are to be included and **false** otherwise.
528 ***********************************************************************************************************/
529
530bool ProtareSingle::useMultiGroupSummedDelayedNeutronsData( Transporting::MG const &a_settings ) const {
531
532 return( ( a_settings.delayedNeutrons( ) == Transporting::DelayedNeutrons::on ) && ( m_multiGroupSummedDelayedNeutrons != nullptr ) );
533}
534
535/* *********************************************************************************************************//**
536 * Returns *a_reaction* except when Rutherford scattering present and only nuclear + Coulomb interference is wanted. Otherwise
537 * returns *m_nuclearPlusCoulombInterferenceOnlyReaction* which may be a **nullptr**.
538 *
539 * @param a_settings [in] Specifies user options.
540 * @param a_reaction [in] Reaction pointer to return if check passes.
541 *
542 * @return Const reaction pointer or *m_nuclearPlusCoulombInterferenceOnlyReaction*.
543 ***********************************************************************************************************/
544
546
547 if( a_reaction->RutherfordScatteringPresent( ) && a_settings.nuclearPlusCoulombInterferenceOnly( ) ) {
548 return( m_nuclearPlusCoulombInterferenceOnlyReaction );
549 }
550
551 return( a_reaction );
552}
553
554/* *********************************************************************************************************//**
555 * Returns *a_reaction* except when Rutherford scattering present and only nuclear + Coulomb interference is wanted. Otherwise
556 * returns *m_nuclearPlusCoulombInterferenceOnlyReaction* which may be a **nullptr**.
557 *
558 * @param a_settings [in] Specifies user options.
559 * @param a_reaction [in] Reaction pointer to return if check passes.
560 *
561 * @return Const reaction pointer or *m_nuclearPlusCoulombInterferenceOnlyReaction*.
562 ***********************************************************************************************************/
563
564Reaction const *ProtareSingle::reactionToMultiGroup( Transporting::MG const &a_settings, std::size_t a_index,
565 ExcludeReactionsSet const &a_reactionsToExclude ) const {
566
567 Reaction const *reaction1 = m_reactions.get<Reaction>( a_index );
568
569 if( !reaction1->active( ) ) return( nullptr );
570 if( a_reactionsToExclude.find( a_index ) != a_reactionsToExclude.end( ) ) return( nullptr );
571
572 return( checkIf_nuclearPlusCoulombInterferenceWanted( a_settings, reaction1 ) );
573}
574
575/* *********************************************************************************************************//**
576 * Returns the pointer representing the protare (i.e., *this*) if *a_index* is 0 and nullptr otherwise.
577 *
578 * @param a_index [in] Must always be 0.
579 *
580 * @return Returns the pointer representing *this*.
581 ***********************************************************************************************************/
582
583ProtareSingle *ProtareSingle::protare( std::size_t a_index ) {
584
585 if( a_index != 0 ) return( nullptr );
586 return( this );
587}
588
589/* *********************************************************************************************************//**
590 * Returns the pointer representing the protare (i.e., *this*) if *a_index* is 0 and nullptr otherwise.
591 *
592 * @param a_index [in] Must always be 0.
593 *
594 * @return Returns the pointer representing *this*.
595 ***********************************************************************************************************/
596
597ProtareSingle const *ProtareSingle::protare( std::size_t a_index ) const {
598
599 if( a_index != 0 ) return( nullptr );
600 return( this );
601}
602
603/* *********************************************************************************************************//**
604 * Returns the intid for the requested particle or -1 if the particle is not in *this* PoPs database.
605 *
606 * @param a_id [in] The GNDS PoPs id for particle whose intd is requested.
607 *
608 * @return C++ int for the requested particle or -1 if particle is not in PoPs.
609 ******************************************************************/
610
611int ProtareSingle::intid( std::string const &a_id ) const {
612
613 return( m_internalPoPs.intid( a_id ) );
614}
615
616/* *********************************************************************************************************//**
617 * Fills in a std::set with a unique list of all product indices produced by reactions and orphanProducts.
618 * If a_transportablesOnly is true, only transportable product incides are return.
619 *
620 * @param a_ids [out] The unique list of product indices.
621 * @param a_particles [in] The list of particles to be transported.
622 * @param a_transportablesOnly [in] If true, only transportable product indices are added in the list.
623 ***********************************************************************************************************/
624
625void ProtareSingle::productIDs( std::set<std::string> &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly ) const {
626
627 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
628 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
629
630 if( !reaction1->active( ) ) continue;
631 reaction1->productIDs( a_ids, a_particles, a_transportablesOnly );
632 }
633
634 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
635 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
636
637 if( !reaction1->active( ) ) continue;
638 reaction1->productIDs( a_ids, a_particles, a_transportablesOnly );
639 }
640}
641
642/* *********************************************************************************************************//**
643 * Determines the maximum Legredre order present in the multi-group transfer matrix for a give product for a give label.
644 *
645 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
646 * @param a_settings [in] Specifies the requested label.
647 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
648 * @param a_productID [in] The id of the requested product.
649 *
650 * @return The maximum Legredre order. If no transfer matrix data are present for the requested product, -1 is returned.
651 ***********************************************************************************************************/
652
654 Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID ) const {
655
656 int _maximumLegendreOrder = -1;
657 ExcludeReactionsSet excludeReactionsSet;
658
659 if( useMultiGroupSummedData( a_settings, excludeReactionsSet ) ) {
660 _maximumLegendreOrder = m_multiGroupSummedReaction->maximumLegendreOrder( a_smr, a_settings, a_temperatureInfo, a_productID );
661 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
662 _maximumLegendreOrder = std::max( _maximumLegendreOrder, m_multiGroupSummedDelayedNeutrons->maximumLegendreOrder(
663 a_smr, a_settings, a_temperatureInfo, a_productID ) );
664 } }
665 else {
666 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
667 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
668
669 if( !reaction1->active( ) ) continue;
670 int r_maximumLegendreOrder = reaction1->maximumLegendreOrder( a_smr, a_settings, a_temperatureInfo, a_productID );
671 if( r_maximumLegendreOrder > _maximumLegendreOrder ) _maximumLegendreOrder = r_maximumLegendreOrder;
672 }
673
674 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
675 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
676
677 if( !reaction1->active( ) ) continue;
678 int r_maximumLegendreOrder = reaction1->maximumLegendreOrder( a_smr, a_settings, a_temperatureInfo, a_productID );
679 if( r_maximumLegendreOrder > _maximumLegendreOrder ) _maximumLegendreOrder = r_maximumLegendreOrder;
680 }
681 }
682
683 return( _maximumLegendreOrder );
684}
685
686/* *********************************************************************************************************//**
687 * Returns a list of all process temperature data. For each temeprature, the labels for its
688 *
689 * + heated cross section data,
690 * + gridded cross section data,
691 * + multi-group data, and
692 * + multi-group upscatter data.
693 *
694 * are returned. If no data are present for a give data type (e.g., gridded cross section, multi-group upscatter), its label is an empty std::string.
695 *
696 * @return The list of temperatures and their labels via an Styles::TemperatureInfos instance. The Styles::TemperatureInfos class
697 * has many (if not all) the method of a std::vector.
698 ***********************************************************************************************************/
699
701
702 std::size_t size( m_styles.size( ) );
703 Styles::TemperatureInfos temperature_infos;
704
705 for( std::size_t i1 = 0; i1 < size; ++i1 ) {
706 Styles::Base const *style1 = m_styles.get<Styles::Base>( i1 );
707
708 if( style1->moniker( ) == GIDI_heatedStyleChars ) {
709 PhysicalQuantity const &temperature = style1->temperature( );
710 std::string heated_cross_section( style1->label( ) );
711 std::string gridded_cross_section( "" );
712 std::string URR_probability_tables( "" );
713 std::string heated_multi_group( "" );
714 std::string Sn_elastic_upscatter( "" );
715
716 for( std::size_t i2 = 0; i2 < size; ++i2 ) {
717 Styles::Base const *style2 = m_styles.get<Styles::Base>( i2 );
718
719 if( style2->moniker( ) == GIDI_multiGroupStyleChars ) continue;
720 if( style2->temperature( ).value( ) != temperature.value( ) ) continue;
721
722 if( style2->moniker( ) == GIDI_griddedCrossSectionStyleChars ) {
723 gridded_cross_section = style2->label( ); }
724 else if( style2->moniker( ) == GIDI_URR_probabilityTablesStyleChars ) {
725 URR_probability_tables = style2->label( ); }
726 else if( style2->moniker( ) == GIDI_SnElasticUpScatterStyleChars ) {
727 Sn_elastic_upscatter = style2->label( ); }
728 else if( style2->moniker( ) == GIDI_heatedMultiGroupStyleChars ) {
729 heated_multi_group = style2->label( );
730 }
731 }
732 temperature_infos.push_back( Styles::TemperatureInfo( temperature, heated_cross_section, gridded_cross_section, URR_probability_tables,
733 heated_multi_group, Sn_elastic_upscatter ) );
734 }
735 }
736
737 std::sort( temperature_infos.begin( ), temperature_infos.end( ), sortTemperatures );
738
739 return( temperature_infos );
740}
741
742/* *********************************************************************************************************//**
743 * FOR INTERNAL USE ONLY.
744 *
745 * Determines if the temperature of lhs is less than that of rhs, or not.
746 *
747 * @param lhs [in]
748 * @param rhs [in]
749 * @return true if temperature of lhs < rhs and false otherwise.
750 ***********************************************************************************************************/
751
752bool sortTemperatures( Styles::TemperatureInfo const &lhs, Styles::TemperatureInfo const &rhs ) {
753
754 return( lhs.temperature( ).value( ) < rhs.temperature( ).value( ) );
755}
756
757/* *********************************************************************************************************//**
758 * Returns the (*a_index*+1)th reaction or **nullptr**. If the indexed reaction is deactivated or exlucded,
759 * a **nullptr** is returned.
760 *
761 * @param a_index [in] The index of the requested reaction.
762 * @param a_settings [in] Specifies the requested label.
763 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the cross section.
764 *
765 * @return The (*a_index*+1)th reaction or **nullptr**.
766 ***********************************************************************************************************/
767
768Reaction const *ProtareSingle::reaction( std::size_t a_index, Transporting::MG const &a_settings,
769 ExcludeReactionsSet const &a_reactionsToExclude ) const {
770
771 return( reactionToMultiGroup( a_settings, a_index, a_reactionsToExclude ) );
772}
773
774/* *********************************************************************************************************//**
775 * The method returns the number of reactions of *this* that have been deactivated.
776 *
777 * @return The number of deactivated reaction of *this*.
778 ***********************************************************************************************************/
779
781
782 std::size_t count = 0;
783
784 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
785 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
786
787 if( !reaction1->active( ) ) ++count;
788 }
789
790 return( count );
791}
792
793/* *********************************************************************************************************//**
794 * Re-indexs the reactions in the reactions, orphanProducts and fissionComponents suites.
795 *
796 ***********************************************************************************************************/
797
798void ProtareSingle::updateReactionIndices( std::size_t a_offset ) const {
799
800 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
801 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 + a_offset );
802
803 reaction1->setReactionIndex( i1 );
804 }
805 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
806 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
807
808 reaction1->setReactionIndex( i1 );
809 }
810 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
811 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
812
813 reaction1->setReactionIndex( i1 );
814 }
815}
816
817/* *********************************************************************************************************//**
818 * Returns the multi-group boundaries for the requested label and product.
819 *
820 * @param a_settings [in] Specifies the requested label.
821 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
822 * @param a_productID [in] ID for the requested product.
823 *
824 * @return List of multi-group boundaries.
825 ***********************************************************************************************************/
826
827std::vector<double> ProtareSingle::groupBoundaries( LUPI_maybeUnused Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID ) const {
828
829 Styles::HeatedMultiGroup const *heatedMultiGroupStyle1 = m_styles.get<Styles::HeatedMultiGroup>( a_temperatureInfo.heatedMultiGroup( ) );
830
831 return( heatedMultiGroupStyle1->groupBoundaries( a_productID ) );
832}
833
834/* *********************************************************************************************************//**
835 * Returns the inverse speeds for the requested label. The label must be for a heated multi-group style.
836 *
837 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
838 * @param a_settings [in] Specifies the requested label.
839 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
840 *
841 * @return List of inverse speeds.
842 ***********************************************************************************************************/
843
845 Styles::TemperatureInfo const &a_temperatureInfo ) const {
846
847 Styles::HeatedMultiGroup const *heatedMultiGroupStyle1 = m_styles.get<Styles::HeatedMultiGroup>( a_temperatureInfo.heatedMultiGroup( ) );
848
849 return( heatedMultiGroupStyle1->inverseSpeedData( ) );
850}
851
852/* *********************************************************************************************************//**
853 * Returns true if at least one reaction contains a fission channel.
854 *
855 * @return true if at least one reaction contains a fission channel.
856 ***********************************************************************************************************/
857
859
860 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
861 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
862
863 if( !reaction1->active( ) ) continue;
864 if( reaction1->hasFission( ) ) return( true );
865 }
866 return( false );
867}
868
869/* *********************************************************************************************************//**
870 * Returns **false* if protare has delayed fission neutrons for an active reaction and they are not complete; otherwise, returns **true**.
871 *
872 * @return bool
873 ***********************************************************************************************************/
874
876
877 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
878 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
879
880 if( !reaction1->active( ) ) continue;
881 if( reaction1->hasFission( ) ) {
882 OutputChannel const *outputChannel = reaction1->outputChannel( );
883 if( !outputChannel->isDelayedFissionNeutronComplete( ) ) return( false );
884 }
885 }
886
887 return( true );
888}
889
890/* *********************************************************************************************************//**
891 * 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.
892 *
893 * @param a_item [in] The name of the class member whose pointer is to be return.
894 * @return The pointer to the class member or nullptr if class does not have a member named a_item.
895 ***********************************************************************************************************/
896
897GUPI::Ancestry *ProtareSingle::findInAncestry3( std::string const &a_item ) {
898
899 if( a_item == GIDI_stylesChars ) return( &m_styles );
900 if( a_item == GIDI_reactionsChars ) return( &m_reactions );
901 if( a_item == GIDI_orphanProductsChars ) return( &m_orphanProducts );
902 if( a_item == GIDI_sumsChars ) return( &m_sums );
903 if( a_item == GIDI_fissionComponentsChars ) return( &m_fissionComponents );
904 if( a_item == GIDI_ACE_URR_probabilityTablesChars ) return( &m_ACE_URR_probabilityTables );
905
906 return( nullptr );
907}
908
909/* *********************************************************************************************************//**
910 * 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.
911 *
912 * @param a_item [in] The name of the class member whose pointer is to be return.
913 * @return The pointer to the class member or nullptr if class does not have a member named a_item.
914 ***********************************************************************************************************/
915
916GUPI::Ancestry const *ProtareSingle::findInAncestry3( std::string const &a_item ) const {
917
918 if( a_item == GIDI_stylesChars ) return( &m_styles );
919 if( a_item == GIDI_reactionsChars ) return( &m_reactions );
920 if( a_item == GIDI_orphanProductsChars ) return( &m_orphanProducts );
921 if( a_item == GIDI_sumsChars ) return( &m_sums );
922 if( a_item == GIDI_fissionComponentsChars ) return( &m_fissionComponents );
923 if( a_item == GIDI_ACE_URR_probabilityTablesChars ) return( &m_ACE_URR_probabilityTables );
924
925 return( nullptr );
926}
927
928/* *********************************************************************************************************//**
929 * Returns the multi-group, total cross section for the requested label. This is summed over all reactions.
930 *
931 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
932 * @param a_settings [in] Specifies the requested label.
933 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
934 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the cross section.
935 * @param 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.
936 *
937 * @return The requested multi-group cross section as a GIDI::Vector.
938 ***********************************************************************************************************/
939
941 Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude, std::string const &a_label ) const {
942
943 Vector vector;
944
945 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
946 vector = m_multiGroupSummedReaction->multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo, a_label ); }
947 else {
948 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
949 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
950
951 if( reaction1 != nullptr ) vector += reaction1->multiGroupCrossSection( a_smr, a_settings, a_temperatureInfo, a_label );
952 }
953 }
954
955 return( vector );
956}
957
958/* *********************************************************************************************************//**
959 * Returns the multi-group, total multiplicity for the requested label for the requested product. This is a cross section weighted multiplicity.
960 *
961 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
962 * @param a_settings [in] Specifies the requested label.
963 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
964 * @param a_productID [in] Id for the requested product.
965 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the multiplicity.
966 *
967 * @return The requested multi-group multiplicity as a GIDI::Vector.
968 ***********************************************************************************************************/
969
971 Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID,
972 ExcludeReactionsSet const &a_reactionsToExclude ) const {
973
974 Vector vector;
975
976 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
977 vector = m_multiGroupSummedReaction->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID );
978 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
979 vector += m_multiGroupSummedDelayedNeutrons->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID );
980 } }
981 else {
982 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
983 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
984
985 if( reaction1 != nullptr ) vector += reaction1->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID );
986 }
987
988 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
989 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
990
991 if( !reaction1->active( ) ) continue;
992 vector += reaction1->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, a_productID );
993 }
994 }
995
996 return( vector );
997}
998
999/* *********************************************************************************************************//**
1000 * Returns the multi-group, total fission neutron multiplicity for the requested label. This is a cross section weighted multiplicity.
1001 *
1002 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1003 * @param a_settings [in] Specifies the requested label.
1004 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1005 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the cross multiplicity.
1006 *
1007 * @return The requested multi-group fission neutron multiplicity as a GIDI::Vector.
1008 ***********************************************************************************************************/
1009
1011 Styles::TemperatureInfo const &a_temperatureInfo, LUPI_maybeUnused ExcludeReactionsSet const &a_reactionsToExclude ) const {
1012
1013 Vector vector( 0 );
1014
1015 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1016 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
1017
1018 if( !reaction1->active( ) ) continue;
1019 if( reaction1->hasFission( ) ) vector += reaction1->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, PoPI::IDs::neutron );
1020 }
1021
1022 return( vector );
1023}
1024
1025/* *********************************************************************************************************//**
1026 * Returns the multi-group, total fission gamma multiplicity for the requested label. This is a cross section weighted multiplicity.
1027 *
1028 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1029 * @param a_settings [in] Specifies the requested label.
1030 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1031 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the multiplicity.
1032 *
1033 * @return The requested multi-group fission neutron multiplicity as a GIDI::Vector.
1034 ***********************************************************************************************************/
1035
1037 Styles::TemperatureInfo const &a_temperatureInfo, LUPI_maybeUnused ExcludeReactionsSet const &a_reactionsToExclude ) const {
1038
1039 Vector vector( 0 );
1040
1041 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1042 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
1043
1044 if( !reaction1->active( ) ) continue;
1045 if( reaction1->hasFission( ) ) vector += reaction1->multiGroupMultiplicity( a_smr, a_settings, a_temperatureInfo, PoPI::IDs::photon );
1046 }
1047
1048 return( vector );
1049}
1050
1051/* *********************************************************************************************************//**
1052 * Returns the multi-group, total Q for the requested label. This is a cross section weighted Q
1053 * summed over all reactions
1054 *
1055 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1056 * @param a_settings [in] Specifies the requested label.
1057 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1058 * @param a_final [in] If false, only the Q for the primary reactions are return, otherwise, the Q for the final reactions.
1059 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the Q.
1060 *
1061 * @return The requested multi-group Q as a GIDI::Vector.
1062 ***********************************************************************************************************/
1063
1065 Styles::TemperatureInfo const &a_temperatureInfo, bool a_final, LUPI_maybeUnused bool a_effectivePhotoAtomic,
1066 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1067
1068 Vector vector( 0 );
1069
1070 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1071 vector = m_multiGroupSummedReaction->multiGroupQ( a_smr, a_settings, a_temperatureInfo, a_final );
1072 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
1073 vector += m_multiGroupSummedDelayedNeutrons->multiGroupQ( a_smr, a_settings, a_temperatureInfo, a_final );
1074 } }
1075 else {
1076 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1077 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1078
1079 if( reaction1 != nullptr ) vector += reaction1->multiGroupQ( a_smr, a_settings, a_temperatureInfo, a_final );
1080 }
1081 }
1082
1083 return( vector );
1084}
1085
1086/* *********************************************************************************************************//**
1087 * Returns the multi-group, total product matrix for the requested label for the requested product id for the requested Legendre order.
1088 * If no data are found, an empty GIDI::Matrix is returned.
1089 *
1090 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1091 * @param a_settings [in] Specifies the requested label and if delayed neutrons should be included.
1092 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1093 * @param a_particles [in] The list of particles to be transported.
1094 * @param a_productID [in] PoPs id for the requested product.
1095 * @param a_order [in] Requested product matrix, Legendre order.
1096 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the product matrix.
1097 *
1098 * @return The requested multi-group product matrix as a GIDI::Matrix.
1099 ***********************************************************************************************************/
1100
1102 Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles,
1103 std::string const &a_productID, std::size_t a_order, ExcludeReactionsSet const &a_reactionsToExclude ) const {
1104
1105 Matrix matrix( 0, 0 );
1106
1107 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1108 matrix = m_multiGroupSummedReaction->multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_productID, a_order );
1109 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
1110 matrix += m_multiGroupSummedDelayedNeutrons->multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_productID, a_order );
1111 } }
1112 else {
1113 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1114 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1115
1116 if( reaction1 != nullptr ) matrix += reaction1->multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_productID, a_order );
1117 }
1118
1119 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
1120 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
1121
1122 if( !reaction1->active( ) ) continue;
1123 matrix += reaction1->multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_productID, a_order );
1124 }
1125 }
1126
1127 return( matrix );
1128}
1129
1130/* *********************************************************************************************************//**
1131 * Like ProtareSingle::multiGroupProductMatrix, but only returns the fission neutron, transfer matrix.
1132 *
1133 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1134 * @param a_settings [in] Specifies the requested label and if delayed neutrons should be included.
1135 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1136 * @param a_particles [in] The list of particles to be transported.
1137 * @param a_order [in] Requested product matrix, Legendre order.
1138 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the fission matrix.
1139 *
1140 * @return The requested multi-group neutron fission matrix as a GIDI::Matrix.
1141 ***********************************************************************************************************/
1142
1144 Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::size_t a_order,
1145 LUPI_maybeUnused ExcludeReactionsSet const &a_reactionsToExclude ) const {
1146
1147 Matrix matrix( 0, 0 );
1148
1149 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1150 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
1151
1152 if( !reaction1->active( ) ) continue;
1153 if( reaction1->hasFission( ) ) matrix += reaction1->multiGroupFissionMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, a_order );
1154 }
1155
1156 return( matrix );
1157}
1158
1159/* *********************************************************************************************************//**
1160 * Returns the multi-group transport correction for the requested label. The transport correction is calculated from the transfer matrix
1161 * for the projectile id for the Legendre order of *a_order + 1*.
1162 *
1163 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1164 * @param a_settings [in] Specifies the requested label.
1165 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1166 * @param a_particles [in] The list of particles to be transported.
1167 * @param a_order [in] Maximum Legendre order for transport. The returned transport correction is for the next higher Legender order.
1168 * @param a_transportCorrectionType [in] Requested transport correction type.
1169 * @param a_temperature [in] The temperature of the flux to use when collapsing. Pass to the GIDI::collapse method.
1170 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the transport correction.
1171 *
1172 * @return The requested multi-group transport correction as a GIDI::Vector.
1173 ***********************************************************************************************************/
1174
1176 Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::size_t a_order,
1177 TransportCorrectionType a_transportCorrectionType, double a_temperature, LUPI_maybeUnused ExcludeReactionsSet const &a_reactionsToExclude ) const {
1178
1179 if( a_transportCorrectionType == TransportCorrectionType::None ) return( Vector( 0 ) );
1180
1181 Matrix matrix( multiGroupProductMatrix( a_smr, a_settings, a_temperatureInfo, a_particles, projectile( ).ID( ), a_order + 1 ) );
1182 Matrix matrixCollapsed = collapse( matrix, a_settings, a_particles, a_temperature, projectile( ).ID( ) );
1183 std::size_t size = matrixCollapsed.size( );
1184 std::vector<double> transportCorrection1( size, 0 );
1185
1186 if( a_transportCorrectionType == TransportCorrectionType::None ) {
1187 }
1188 else if( a_transportCorrectionType == TransportCorrectionType::Pendlebury ) {
1189 for( std::size_t index = 0; index < size; ++index ) transportCorrection1[index] = matrixCollapsed[index][index]; }
1190 else {
1191 throw Exception( "Unsupported transport correction: only None and Pendlebury (i.e., Pendlebury/Underhill) are currently supported." );
1192 }
1193 return( Vector( transportCorrection1 ) );
1194}
1195
1196/* *********************************************************************************************************//**
1197 * Returns the multi-group, total available energy for the requested label. This is a cross section weighted available energy
1198 * summed over all reactions.
1199 *
1200 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1201 * @param a_settings [in] Specifies the requested label.
1202 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1203 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the available energy.
1204 *
1205 * @return The requested multi-group available energy as a GIDI::Vector.
1206 ***********************************************************************************************************/
1207
1209 Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude ) const {
1210
1211 Vector vector( 0 );
1212
1213 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1214 vector = m_multiGroupSummedReaction->multiGroupAvailableEnergy( a_smr, a_settings, a_temperatureInfo ); }
1215 else {
1216 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1217 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1218
1219 if( reaction1 != nullptr ) vector += reaction1->multiGroupAvailableEnergy( a_smr, a_settings, a_temperatureInfo );
1220 }
1221 }
1222
1223 return( vector );
1224}
1225
1226/* *********************************************************************************************************//**
1227 * Returns the multi-group, total average energy for the requested label for the requested product. This is a cross section weighted average energy
1228 * summed over all reactions.
1229 *
1230 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1231 * @param a_settings [in] Specifies the requested label.
1232 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1233 * @param a_productID [in] Particle id for the requested product.
1234 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the average energy.
1235 *
1236 * @return The requested multi-group average energy as a GIDI::Vector.
1237 ***********************************************************************************************************/
1238
1240 Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID,
1241 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1242
1243 Vector vector( 0 );
1244
1245 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1246 vector = m_multiGroupSummedReaction->multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, a_productID );
1247 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
1248 vector += m_multiGroupSummedDelayedNeutrons->multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, a_productID );
1249 } }
1250 else {
1251 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1252 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1253
1254 if( reaction1 != nullptr ) vector += reaction1->multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, a_productID );
1255 }
1256
1257 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
1258 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
1259
1260 if( !reaction1->active( ) ) continue;
1261 vector += reaction1->multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, a_productID );
1262 }
1263 }
1264
1265 return( vector );
1266}
1267
1268/* *********************************************************************************************************//**
1269 * Returns the multi-group, total deposition energy for the requested label. This is a cross section weighted deposition energy
1270 * summed over all reactions. The deposition energy is calculated by subtracting the average energy from each transportable particle
1271 * from the available energy. The list of transportable particles is specified via the list of particle specified in the *a_settings* argument.
1272 *
1273 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1274 * @param a_settings [in] Specifies the requested label and the products that are transported.
1275 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1276 * @param a_particles [in] The list of particles to be transported.
1277 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the deposition energy.
1278 *
1279 * @return The requested multi-group deposition energy as a GIDI::Vector.
1280 ***********************************************************************************************************/
1281
1283 Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles,
1284 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1285
1286 bool atLeastOneReactionHasAllParticlesTracked = false;
1287 Vector vector;
1288
1289 if( a_settings.zeroDepositionIfAllProductsTracked( ) ) {
1290 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1291 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1292
1293 if( reaction1 != nullptr ) {
1294 if( reaction1->areAllProductsTracked( a_particles ) ) {
1295 atLeastOneReactionHasAllParticlesTracked = true;
1296 break;
1297 }
1298 }
1299 }
1300 }
1301
1302 if( atLeastOneReactionHasAllParticlesTracked ) {
1303 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1304 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1305
1306 if( reaction1 != nullptr ) vector += reaction1->multiGroupDepositionEnergy( a_smr, a_settings, a_temperatureInfo, a_particles );
1307 }
1308 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
1309 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
1310
1311 if( !reaction1->active( ) ) continue;
1312 vector += reaction1->multiGroupDepositionEnergy( a_smr, a_settings, a_temperatureInfo, a_particles );
1313 } }
1314 else {
1315 std::map<std::string, Transporting::Particle> const &products( a_particles.particles( ) );
1316 vector = multiGroupAvailableEnergy( a_smr, a_settings, a_temperatureInfo, a_reactionsToExclude );
1317 Vector availableEnergy( vector );
1318
1319 for( std::map<std::string, Transporting::Particle>::const_iterator iter = products.begin( ); iter != products.end( ); ++iter ) {
1320 vector -= multiGroupAverageEnergy( a_smr, a_settings, a_temperatureInfo, iter->first, a_reactionsToExclude );
1321 }
1322
1323 for( std::size_t index = 0; index < availableEnergy.size( ); ++index ) { // Check for values that should probably be 0.0.
1324 if( std::fabs( vector[index] ) < std::fabs( 1e-14 * availableEnergy[index] ) ) vector[index] = 0.0;
1325 }
1326 }
1327
1328 return( vector );
1329}
1330
1331/* *********************************************************************************************************//**
1332 * Returns the multi-group, total available momentum for the requested label. This is a cross section weighted available momentum
1333 * summed over all reactions.
1334 *
1335 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1336 * @param a_settings [in] Specifies the requested label.
1337 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1338 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the available momentum.
1339 *
1340 * @return The requested multi-group available momentum as a GIDI::Vector.
1341 ***********************************************************************************************************/
1342
1344 Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude ) const {
1345
1346 Vector vector( 0 );
1347
1348 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1349 vector = m_multiGroupSummedReaction->multiGroupAvailableMomentum( a_smr, a_settings, a_temperatureInfo ); }
1350 else {
1351 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1352 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1353
1354 if( reaction1 != nullptr ) vector += reaction1->multiGroupAvailableMomentum( a_smr, a_settings, a_temperatureInfo );
1355 }
1356 }
1357
1358 return( vector );
1359}
1360
1361/* *********************************************************************************************************//**
1362 * Returns the multi-group, total average momentum for the requested label for the requested product. This is a cross section weighted average momentum
1363 * summed over all reactions.
1364 *
1365 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1366 * @param a_settings [in] Specifies the requested label.
1367 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1368 * @param a_productID [in] Particle id for the requested product.
1369 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the average momentum.
1370 *
1371 * @return The requested multi-group average momentum as a GIDI::Vector.
1372 ***********************************************************************************************************/
1373
1375 Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID,
1376 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1377
1378 Vector vector( 0 );
1379
1380 if( useMultiGroupSummedData( a_settings, a_reactionsToExclude ) ) {
1381 vector = m_multiGroupSummedReaction->multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, a_productID );
1382 if( useMultiGroupSummedDelayedNeutronsData( a_settings ) ) {
1383 vector += m_multiGroupSummedDelayedNeutrons->multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, a_productID );
1384 } }
1385 else {
1386 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1387 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1388
1389 if( reaction1 != nullptr ) vector += reaction1->multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, a_productID );
1390 }
1391
1392 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
1393 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
1394
1395 if( !reaction1->active( ) ) continue;
1396 vector += reaction1->multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, a_productID );
1397 }
1398 }
1399
1400 return( vector );
1401}
1402
1403/* *********************************************************************************************************//**
1404 * Returns the multi-group, total deposition momentum for the requested label. This is a cross section weighted deposition momentum
1405 * summed over all reactions. The deposition momentum is calculated by subtracting the average momentum from each transportable particle
1406 * from the available momentum. The list of transportable particles is specified via the list of particle specified in the *a_settings* argument.
1407 *
1408 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1409 * @param a_settings [in] Specifies the requested label.
1410 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1411 * @param a_particles [in] The list of particles to be transported.
1412 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the deposition momentum.
1413 *
1414 * @return The requested multi-group deposition momentum as a GIDI::Vector.
1415 ***********************************************************************************************************/
1416
1418 Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles,
1419 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1420
1421 std::map<std::string, Transporting::Particle> const &products( a_particles.particles( ) );
1422 Vector vector = multiGroupAvailableMomentum( a_smr, a_settings, a_temperatureInfo, a_reactionsToExclude );
1423
1424 for( std::map<std::string, Transporting::Particle>::const_iterator iter = products.begin( ); iter != products.end( ); ++iter ) {
1425 vector -= multiGroupAverageMomentum( a_smr, a_settings, a_temperatureInfo, iter->first, a_reactionsToExclude );
1426 }
1427
1428 return( vector );
1429}
1430
1431/* *********************************************************************************************************//**
1432 * Returns the multi-group, gain for the requested particle and label. This is a cross section weighted gain summed over all reactions.
1433 *
1434 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1435 * @param a_settings [in] Specifies the requested label.
1436 * @param a_temperatureInfo [in] Specifies the temperature and labels use to lookup the requested data.
1437 * @param a_productID [in] The particle PoPs' id for the whose gain is to be calculated.
1438 * @param a_reactionsToExclude [in] A list of reaction indices that are to be ignored when calculating the gain.
1439 *
1440 * @return The requested multi-group gain as a **GIDI::Vector**.
1441 ***********************************************************************************************************/
1442
1444 Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID,
1445 ExcludeReactionsSet const &a_reactionsToExclude ) const {
1446
1447 Vector vector( 0 );
1448 std::string const projectile_ID = projectile( ).ID( );
1449
1450 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1451 Reaction const *reaction1 = reactionToMultiGroup( a_settings, i1, a_reactionsToExclude );
1452
1453 if( reaction1 != nullptr ) vector += reaction1->multiGroupGain( a_smr, a_settings, a_temperatureInfo, a_productID, projectile_ID );
1454 }
1455
1456 return( vector );
1457}
1458
1459/* *********************************************************************************************************//**
1460 *
1461 *
1462 * @return A list of label, mu cutoff pairs.
1463 ***********************************************************************************************************/
1464
1466
1467 stringAndDoublePairs muCutoffs;
1468
1469 for( std::size_t i1 = 0; i1 < m_styles.size( ); ++i1 ) {
1470 Styles::Base const *style1 = m_styles.get<Styles::Base>( i1 );
1471
1473 Styles::CoulombPlusNuclearElasticMuCutoff const *style2 = static_cast<Styles::CoulombPlusNuclearElasticMuCutoff const *>( style1 );
1474
1475 stringAndDoublePair labelMu( style2->label( ), style2->muCutoff( ) );
1476
1477 muCutoffs.push_back( std::move( labelMu ) );
1478 }
1479 }
1480
1481 return( muCutoffs );
1482}
1483
1484/* *********************************************************************************************************//**
1485 * Returns the list of DelayedNeutronProduct instances.
1486 *
1487 * @return a_delayedNeutronProducts The list of delayed neutrons.
1488 ***********************************************************************************************************/
1489
1491
1492 DelayedNeutronProducts delayedNeutronProducts1;
1493
1494 if( hasFission( ) ) {
1495 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1496 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
1497
1498 if( !reaction1->active( ) ) continue;
1499 if( reaction1->hasFission( ) ) reaction1->delayedNeutronProducts( delayedNeutronProducts1 );
1500 }
1501 }
1502
1503 return( delayedNeutronProducts1 );
1504}
1505
1506/* *********************************************************************************************************//**
1507 * Calls the **incompleteParticles** method for each active reaction in the *reactions* and *orphanProducts* nodes.
1508 *
1509 * @param a_settings [in] Specifies the requested label.
1510 * @param a_incompleteParticles [out] The list of particles whose **completeParticle** method returns *false*.
1511 ***********************************************************************************************************/
1512
1513void ProtareSingle::incompleteParticles( Transporting::Settings const &a_settings, std::set<std::string> &a_incompleteParticles ) const {
1514
1515 for( std::size_t i1 = 0; i1 < m_reactions.size( ); ++i1 ) {
1516 Reaction const *reaction1 = m_reactions.get<Reaction>( i1 );
1517
1518 if( !reaction1->active( ) ) continue;
1519 reaction1->incompleteParticles( a_settings, a_incompleteParticles );
1520 }
1521
1522 for( std::size_t i1 = 0; i1 < m_orphanProducts.size( ); ++i1 ) {
1523 Reaction const *reaction1 = m_orphanProducts.get<Reaction>( i1 );
1524
1525 if( !reaction1->active( ) ) continue;
1526 reaction1->incompleteParticles( a_settings, a_incompleteParticles );
1527 }
1528}
1529
1530/* *********************************************************************************************************//**
1531 * Write *this* to a file in GNDS/XML format.
1532 *
1533 * @param a_fileName [in] Name of file to save XML lines to.
1534 ***********************************************************************************************************/
1535
1536void ProtareSingle::saveAs( std::string const &a_fileName ) const {
1537
1538 GUPI::WriteInfo writeInfo;
1539
1540 toXMLList( writeInfo, "" );
1541
1542 std::ofstream fileio;
1543 fileio.open( a_fileName.c_str( ) );
1544 for( std::list<std::string>::iterator iter = writeInfo.m_lines.begin( ); iter != writeInfo.m_lines.end( ); ++iter ) {
1545 fileio << *iter << std::endl;
1546 }
1547 fileio.close( );
1548}
1549
1550/* *********************************************************************************************************//**
1551 * Fills the argument *a_writeInfo* with the XML lines that represent *this*. Recursively enters each sub-node.
1552 *
1553 * @param a_writeInfo [in/out] Instance containing incremental indentation and other information and stores the appended lines.
1554 * @param a_indent [in] The amount to indent *this* node.
1555 ***********************************************************************************************************/
1556
1557void ProtareSingle::toXMLList( GUPI::WriteInfo &a_writeInfo, std::string const &a_indent ) const {
1558
1559 std::string indent2 = a_writeInfo.incrementalIndent( a_indent );
1560 std::string header = LUPI_XML_verionEncoding;
1561 std::string attributes;
1562
1563 a_writeInfo.push_back( header );
1564
1565 attributes = a_writeInfo.addAttribute( GIDI_projectileChars, projectile( ).ID( ) );
1566 attributes += a_writeInfo.addAttribute( GIDI_targetChars, GNDS_target( ).ID( ) );
1567 attributes += a_writeInfo.addAttribute( GIDI_evaluationChars, evaluation( ) );
1568 attributes += a_writeInfo.addAttribute( GIDI_formatChars, m_formatVersion.format( ) );
1569 attributes += a_writeInfo.addAttribute( GIDI_projectileFrameChars, frameToString( projectileFrame( ) ) );
1570
1571 a_writeInfo.addNodeStarter( a_indent, moniker( ), attributes );
1572
1573 m_externalFiles.toXMLList( a_writeInfo, indent2 );
1574 m_styles.toXMLList( a_writeInfo, indent2 );
1575
1576 std::vector<std::string> pops_XMLList;
1577 m_internalPoPs.toXMLList( pops_XMLList, indent2 );
1578 for( std::vector<std::string>::iterator iter = pops_XMLList.begin( ); iter != pops_XMLList.end( ); ++iter ) a_writeInfo.push_back( *iter );
1579
1580 m_reactions.toXMLList( a_writeInfo, indent2 );
1581 m_orphanProducts.toXMLList( a_writeInfo, indent2 );
1582 m_sums.toXMLList( a_writeInfo, indent2 );
1583 m_fissionComponents.toXMLList( a_writeInfo, indent2 );
1584
1585 a_writeInfo.addNodeEnder( moniker( ) );
1586}
1587
1588/* *********************************************************************************************************//**
1589 * This method parses the **targetInfo** node in the **evaluated** style node into the *m_targetInfo* member of *this*.
1590 *
1591 * @param a_node [in] The protare (i.e., reactionSuite) node to be parsed and used to construct a Protare.
1592 ***********************************************************************************************************/
1593
1595
1596 m_targetInfo.parseEvaluatedTargetInfo( a_node );
1597}
1598
1599}
#define GIDI_multiGroupStyleChars
Definition GIDI.hpp:246
#define GIDI_MapInteractionTNSLChars
Definition GIDI.hpp:5165
#define GIDI_MapInteractionAtomicChars
Definition GIDI.hpp:5164
#define GIDI_ENDF_MT_Chars
Definition GIDI.hpp:420
#define GIDI_stylesChars
Definition GIDI.hpp:177
#define GIDI_griddedCrossSectionStyleChars
Definition GIDI.hpp:251
#define GIDI_LLNL_photoAtomicIncoherentDoppler_Chars
Definition GIDI.hpp:188
#define GIDI_CoulombPlusNuclearElasticMuCutoffStyleChars
Definition GIDI.hpp:243
#define GIDI_URR_probabilityTablesStyleChars
Definition GIDI.hpp:252
#define GIDI_projectileChars
Definition GIDI.hpp:413
#define GIDI_outputChannelChars
Definition GIDI.hpp:227
#define GIDI_externalFilesChars
Definition GIDI.hpp:172
#define GIDI_nuclearPlusCoulombInterferenceChars
Definition GIDI.hpp:202
#define GIDI_applicationDataChars
Definition GIDI.hpp:200
#define GIDI_reactionChars
Definition GIDI.hpp:180
#define GIDI_reactionsChars
Definition GIDI.hpp:179
#define GIDI_LLNL_multiGroupReactions_Chars
Definition GIDI.hpp:156
#define GIDI_evaluationChars
Definition GIDI.hpp:415
#define GIDI_LLNL_multiGroupDelayedNeutrons_Chars
Definition GIDI.hpp:157
#define GIDI_PoPsChars
Definition GIDI.hpp:178
#define GIDI_labelChars
Definition GIDI.hpp:438
#define GIDI_LLNL_Chars
Definition GIDI.hpp:155
#define GIDI_heatedMultiGroupStyleChars
Definition GIDI.hpp:253
#define GIDI_heatedStyleChars
Definition GIDI.hpp:250
#define GIDI_sumsChars
Definition GIDI.hpp:204
#define GIDI_orphanProductsChars
Definition GIDI.hpp:181
#define GIDI_doubleDifferentialCrossSectionChars
Definition GIDI.hpp:214
#define GIDI_targetChars
Definition GIDI.hpp:414
#define GIDI_LLNL_URR_probability_tables_Chars
Definition GIDI.hpp:158
#define GIDI_LLNL_GRIN_continuumGammas
Definition GIDI.hpp:160
#define GIDI_GRIN_continuumGammasChars
Definition GIDI.hpp:474
#define GIDI_TNSLChars
Definition GIDI.hpp:165
#define GIDI_LLNL_pointwiseAverageProductEnergies
Definition GIDI.hpp:159
#define GIDI_interactionChars
Definition GIDI.hpp:416
#define GIDI_projectileFrameChars
Definition GIDI.hpp:419
#define GIDI_averageEnergyChars
Definition GIDI.hpp:225
#define GIDI_topLevelChars
Definition GIDI.hpp:169
#define GIDI_SnElasticUpScatterStyleChars
Definition GIDI.hpp:254
#define GIDI_incompleteReactionsChars
Definition GIDI.hpp:183
#define GIDI_formatChars
Definition GIDI.hpp:167
#define GIDI_institutionChars
Definition GIDI.hpp:201
#define GIDI_documentations_1_10_Chars
Definition GIDI.hpp:176
#define GIDI_fissionComponentsChars
Definition GIDI.hpp:184
#define GIDI_ACE_URR_probabilityTablesChars
Definition GIDI.hpp:186
#define LUPI_XML_verionEncoding
Definition LUPI.hpp:28
#define LUPI_maybeUnused
ParseMode parseMode() const
Definition GIDI.hpp:557
std::size_t size() const
bool isDelayedFissionNeutronComplete() const
PhysicalQuantity const & mass() const
Definition GIDI.hpp:765
std::string const & ID() const
Definition GIDI.hpp:760
double value() const
Definition GIDI.hpp:728
Vector multiGroupFissionGammaMultiplicity(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Reaction * reaction(std::size_t a_index)
Definition GIDI.hpp:4816
Styles::TemperatureInfos temperatures() const
void updateReactionIndices(std::size_t a_offset) 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, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupTransportCorrection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, std::size_t a_order, TransportCorrectionType a_transportCorrectionType, double a_temperature, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
std::vector< double > groupBoundaries(Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Vector multiGroupAvailableMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupAvailableEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
void saveAs(std::string const &a_fileName) const
void parseEvaluatedTargetInfo(HAPI::Node const &a_node)
void toXMLList(GUPI::WriteInfo &a_writeInfo, std::string const &a_indent="") const
void productIDs(std::set< std::string > &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly) const
void incompleteParticles(Transporting::Settings const &a_settings, std::set< std::string > &a_incompleteParticles) const
int maximumLegendreOrder(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
std::string const & evaluation(LUPI_maybeUnused std::size_t a_index=0) const
Definition GIDI.hpp:4786
Vector multiGroupAverageMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
bool hasFission() const
PoPI::Database const & internalPoPs() const
Definition GIDI.hpp:4805
Vector multiGroupDepositionEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
GUPI::Ancestry * findInAncestry3(std::string const &a_item)
bool isDelayedFissionNeutronComplete() const
int intid(std::string const &a_id) const
ProtareSingle(PoPI::Database const &a_pops, std::string const &a_projectileID, std::string const &a_targetID, std::string const &a_evaluation, std::string const &a_interaction, std::string const &a_formatVersion=GNDS_formatVersion_1_10Chars)
Vector multiGroupDepositionMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, Transporting::Particles const &a_particles, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
std::size_t numberOfInactiveReactions() const
Reaction const * checkIf_nuclearPlusCoulombInterferenceWanted(Transporting::MG const &a_settings, Reaction const *a_reaction) const
Frame projectileFrame(LUPI_maybeUnused std::size_t a_index=0) const
Definition GIDI.hpp:4788
DelayedNeutronProducts delayedNeutronProducts() const
Vector multiGroupAverageEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupMultiplicity(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
ProtareSingle * protare(std::size_t a_index)
Vector multiGroupQ(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, bool a_final, bool a_effectivePhotoAtomic=true, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}, std::string const &a_label="") const
stringAndDoublePairs muCutoffForCoulombPlusNuclearElastic() const
Vector multiGroupInverseSpeed(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
Reaction const * reactionToMultiGroup(Transporting::MG const &a_settings, std::size_t a_index, ExcludeReactionsSet const &a_reactionsToExclude) 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, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupGain(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
Vector multiGroupFissionNeutronMultiplicity(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}) const
void setTarget(ParticleInfo const &a_target)
Definition GIDI.hpp:4544
virtual std::size_t numberOfReactions() const =0
virtual void TNSL_crossSectionSumCorrection(std::string const &a_label, Functions::XYs1d &a_crossSectionSum)
void setProjectile(ParticleInfo const &a_projectile)
Definition GIDI.hpp:4542
void initialize(HAPI::Node const &a_node, SetupInfo &a_setupInfo, PoPI::Database const &a_pops, PoPI::Database const &a_internalPoPs, bool a_targetRequiredInGlobalPoPs, bool a_requiredInPoPs=true)
ParticleInfo const & target() const
Definition GIDI.hpp:4543
ExcludeReactionsSet reactionIndicesMatchingENDLCValues(std::set< int > const &a_CValues, bool a_checkActiveState=true)
ParticleInfo const & GNDS_target() const
Definition GIDI.hpp:4547
virtual Reaction * reaction(std::size_t a_index)=0
ParticleInfo const & projectile() const
Definition GIDI.hpp:4541
bool areAllProductsTracked(Transporting::Particles const &a_particles) 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
int ENDL_C() const
Definition GIDI.hpp:4285
Vector multiGroupQ(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, bool a_final) const
Vector multiGroupAvailableEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
bool active() const
Definition GIDI.hpp:4280
Vector multiGroupAvailableMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo) const
OutputChannel * outputChannel() const
Definition GIDI.hpp:4307
bool RutherfordScatteringPresent() const
Definition GIDI.hpp:4290
int maximumLegendreOrder(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_label="") 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
bool hasFission() const
void incompleteParticles(Transporting::Settings const &a_settings, std::set< std::string > &a_incompleteParticles) const
Vector multiGroupAverageEnergy(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
void productIDs(std::set< std::string > &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly) const
Vector multiGroupDepositionEnergy(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
Vector multiGroupAverageMomentum(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
Vector multiGroupMultiplicity(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_productID) const
ParticleSubstitution * m_particleSubstitution
Definition GIDI.hpp:598
virtual PhysicalQuantity const & temperature() const =0
std::string const & label() const
Definition GIDI.hpp:3130
Vector inverseSpeedData() const
Definition GIDI.hpp:3329
std::vector< double > groupBoundaries(std::string const &a_ID) const
PhysicalQuantity const & temperature() const
Definition GIDI.hpp:3430
std::string const & heatedMultiGroup() const
Definition GIDI.hpp:3434
bool useMultiGroupSummedData() const
Definition GIDI.hpp:3764
std::map< std::string, Particle > & particles()
Definition GIDI.hpp:3692
bool zeroDepositionIfAllProductsTracked() const
Definition GIDI.hpp:3734
bool nuclearPlusCoulombInterferenceOnly() const
Definition GIDI.hpp:3730
std::size_t size() const
Definition GIDI_data.hpp:79
void setMoniker(std::string const &a_moniker)
Definition GUPI.hpp:103
void setAncestor(Ancestry *a_ancestor)
Definition GUPI.hpp:106
std::string const & moniker() const
Definition GUPI.hpp:102
Ancestry(std::string const &a_moniker, std::string const &a_attribute="")
void push_back(std::string const &a_line)
Definition GUPI.hpp:53
std::list< std::string > m_lines
Definition GUPI.hpp:45
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 name() const
Definition HAPI_Node.cc:136
bool empty() const
Definition HAPI_Node.cc:150
std::string attribute_as_string(const char *a_name) const
Definition HAPI.hpp:179
void to_next_sibling() const
Definition HAPI_Node.cc:112
Node first_child() const
Definition HAPI_Node.cc:82
Node child(const char *name) const
Definition HAPI_Node.cc:72
int attribute_as_int(const char *a_name) const
Definition HAPI.hpp:185
bool exists(std::string const &a_id) const
const char * name(G4int ptype)
std::vector< Styles::TemperatureInfo > TemperatureInfos
Definition GIDI.hpp:3440
Definition GIDI.hpp:32
std::map< std::string, ParticleInfo > ParticleSubstitution
Definition GIDI.hpp:487
std::pair< std::string, double > stringAndDoublePair
Definition GIDI.hpp:485
Form * parseFissionComponent(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)
FileType
Definition GIDI.hpp:148
Form * parseACE_URR_probabilityTables(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)
Frame
Definition GIDI.hpp:146
Form * parseAverageEnergySuite(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)
std::vector< stringAndDoublePair > stringAndDoublePairs
Definition GIDI.hpp:486
std::string frameToString(Frame a_frame)
Definition GIDI_misc.cc:375
Form * parseExternalFilesSuite(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)
Form * parseReaction(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)
Vector collapse(Vector const &a_vector, Transporting::Settings const &a_settings, Transporting::Particles const &a_particles, double a_temperature)
Form * parseOrphanProduct(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 * parseStylesSuite(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)
std::set< std::size_t > ExcludeReactionsSet
Definition GIDI.hpp:47
Frame parseFrame(HAPI::Node const &a_node, SetupInfo &a_setupInfo, std::string const &a_name)
std::vector< DelayedNeutronProduct > DelayedNeutronProducts
Definition GIDI.hpp:4031
TransportCorrectionType
Definition GIDI.hpp:147
Definition GUPI.hpp:20
Definition LUPI.hpp:40
static std::string const photon
Definition PoPI.hpp:162
static std::string const neutron
Definition PoPI.hpp:164