2025-07-01 07:24:18 -04:00
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
< html xmlns = "http://www.w3.org/1999/xhtml" lang = "en-US" >
< head >
< meta http-equiv = "Content-Type" content = "text/xhtml;charset=UTF-8" / >
< meta http-equiv = "X-UA-Compatible" content = "IE=11" / >
< meta name = "generator" content = "Doxygen 1.13.2" / >
< meta name = "viewport" content = "width=device-width, initial-scale=1" / >
< title > GridFire: src/network/include/gridfire/engine/engine_abstract.h Source File< / title >
< link href = "tabs.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "jquery.js" > < / script >
< script type = "text/javascript" src = "dynsections.js" > < / script >
< script type = "text/javascript" src = "clipboard.js" > < / script >
< link href = "navtree.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "navtreedata.js" > < / script >
< script type = "text/javascript" src = "navtree.js" > < / script >
< script type = "text/javascript" src = "resize.js" > < / script >
< script type = "text/javascript" src = "cookie.js" > < / script >
< link href = "search/search.css" rel = "stylesheet" type = "text/css" / >
< script type = "text/javascript" src = "search/searchdata.js" > < / script >
< script type = "text/javascript" src = "search/search.js" > < / script >
< link href = "doxygen.css" rel = "stylesheet" type = "text/css" / >
< link href = "doxygen-awesome.css" rel = "stylesheet" type = "text/css" / >
< link href = "doxygen-awesome-sidebar-only.css" rel = "stylesheet" type = "text/css" / >
< / head >
< body >
< div id = "top" > <!-- do not remove this div, it is closed by doxygen! -->
< div id = "titlearea" >
< table cellspacing = "0" cellpadding = "0" >
< tbody >
< tr id = "projectrow" >
< td id = "projectalign" >
< div id = "projectname" > GridFire< span id = "projectnumber" >   0.0.1a< / span >
< / div >
< div id = "projectbrief" > General Purpose Nuclear Network< / div >
< / td >
< / tr >
< / tbody >
< / table >
< / div >
<!-- end header part -->
<!-- Generated by Doxygen 1.13.2 -->
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699& dn=expat.txt MIT */
var searchBox = new SearchBox("searchBox", "search/",'.html');
/* @license-end */
< / script >
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699& dn=expat.txt MIT */
$(function() { codefold.init(0); });
/* @license-end */
< / script >
< script type = "text/javascript" src = "menudata.js" > < / script >
< script type = "text/javascript" src = "menu.js" > < / script >
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699& dn=expat.txt MIT */
$(function() {
initMenu('',true,false,'search.php','Search',true);
$(function() { init_search(); });
});
/* @license-end */
< / script >
< div id = "main-nav" > < / div >
< / div > <!-- top -->
< div id = "side-nav" class = "ui-resizable side-nav-resizable" >
< div id = "nav-tree" >
< div id = "nav-tree-contents" >
< div id = "nav-sync" class = "sync" > < / div >
< / div >
< / div >
< div id = "splitbar" style = "-moz-user-select:none;"
class="ui-resizable-handle">
< / div >
< / div >
< script type = "text/javascript" >
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699& dn=expat.txt MIT */
$(function(){initNavTree('engine__abstract_8h_source.html',''); initResizable(true); });
/* @license-end */
< / script >
< div id = "doc-content" >
<!-- window showing the filter options -->
< div id = "MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
< / div >
<!-- iframe showing the search results (closed by default) -->
< div id = "MSearchResultsWindow" >
< div id = "MSearchResults" >
< div class = "SRPage" >
< div id = "SRIndex" >
< div id = "SRResults" > < / div >
< div class = "SRStatus" id = "Loading" > Loading...< / div >
< div class = "SRStatus" id = "Searching" > Searching...< / div >
< div class = "SRStatus" id = "NoMatches" > No Matches< / div >
< / div >
< / div >
< / div >
< / div >
< div class = "header" >
< div class = "headertitle" > < div class = "title" > engine_abstract.h< / div > < / div >
< / div > <!-- header -->
< div class = "contents" >
< a href = "engine__abstract_8h.html" > Go to the documentation of this file.< / a > < div class = "fragment" > < div class = "line" > < a id = "l00001" name = "l00001" > < / a > < span class = "lineno" > 1< / span > < span class = "preprocessor" > #pragma once< / span > < / div >
< div class = "line" > < a id = "l00002" name = "l00002" > < / a > < span class = "lineno" > 2< / span > < / div >
< div class = "line" > < a id = "l00003" name = "l00003" > < / a > < span class = "lineno" > 3< / span > < span class = "preprocessor" > #include " < a class = "code" href = "reaction_8h.html" > gridfire/reaction/reaction.h< / a > " < / span > < / div >
2025-07-01 11:40:51 -04:00
< div class = "line" > < a id = "l00004" name = "l00004" > < / a > < span class = "lineno" > 4< / span > < span class = "preprocessor" > #include " < a class = "code" href = "network_8h.html" > gridfire/network.h< / a > " < / span > < / div >
< div class = "line" > < a id = "l00005" name = "l00005" > < / a > < span class = "lineno" > 5< / span > < span class = "preprocessor" > #include " < a class = "code" href = "screening__abstract_8h.html" > gridfire/screening/screening_abstract.h< / a > " < / span > < / div >
< div class = "line" > < a id = "l00006" name = "l00006" > < / a > < span class = "lineno" > 6< / span > < span class = "preprocessor" > #include " < a class = "code" href = "screening__types_8h.html" > gridfire/screening/screening_types.h< / a > " < / span > < / div >
< div class = "line" > < a id = "l00007" name = "l00007" > < / a > < span class = "lineno" > 7< / span > < / div >
< div class = "line" > < a id = "l00008" name = "l00008" > < / a > < span class = "lineno" > 8< / span > < span class = "preprocessor" > #include < vector> < / span > < / div >
< div class = "line" > < a id = "l00009" name = "l00009" > < / a > < span class = "lineno" > 9< / span > < span class = "preprocessor" > #include < unordered_map> < / span > < / div >
< div class = "line" > < a id = "l00010" name = "l00010" > < / a > < span class = "lineno" > 10< / span > < / div >
< div class = "line" > < a id = "l00023" name = "l00023" > < / a > < span class = "lineno" > 23< / span > < / div >
< div class = "foldopen" id = "foldopen00024" data-start = "{" data-end = "}" >
< div class = "line" > < a id = "l00024" name = "l00024" > < / a > < span class = "lineno" > < a class = "line" href = "namespacegridfire.html" > 24< / a > < / span > < span class = "keyword" > namespace < / span > < a class = "code hl_namespace" href = "namespacegridfire.html" > gridfire< / a > {< / div >
< div class = "line" > < a id = "l00025" name = "l00025" > < / a > < span class = "lineno" > 25< / span > < / div >
< div class = "line" > < a id = "l00032" name = "l00032" > < / a > < span class = "lineno" > 32< / span > < span class = "keyword" > template< / span > < < span class = "keyword" > typename< / span > T> < / div >
< div class = "line" > < a id = "l00033" name = "l00033" > < / a > < span class = "lineno" > < a class = "line" href = "conceptgridfire_1_1_is_arithmetic_or_a_d.html" > 33< / a > < / span > < span class = "keyword" > concept < / span > < a class = "code hl_concept" href = "conceptgridfire_1_1_is_arithmetic_or_a_d.html" > IsArithmeticOrAD< / a > = std::is_same_v< T, double> || std::is_same_v< T, CppAD::AD< double> > ;< / div >
< div class = "line" > < a id = "l00034" name = "l00034" > < / a > < span class = "lineno" > 34< / span > < / div >
< div class = "line" > < a id = "l00052" name = "l00052" > < / a > < span class = "lineno" > 52< / span > < span class = "keyword" > template< / span > < IsArithmeticOrAD T> < / div >
< div class = "foldopen" id = "foldopen00053" data-start = "{" data-end = "};" >
< div class = "line" > < a id = "l00053" name = "l00053" > < / a > < span class = "lineno" > < a class = "line" href = "structgridfire_1_1_step_derivatives.html" > 53< / a > < / span > < span class = "keyword" > struct < / span > < a class = "code hl_struct" href = "structgridfire_1_1_step_derivatives.html" > StepDerivatives< / a > {< / div >
< div class = "line" > < a id = "l00054" name = "l00054" > < / a > < span class = "lineno" > < a class = "line" href = "structgridfire_1_1_step_derivatives.html#ae0de268b86c2404379409c4feae0b34d" > 54< / a > < / span > std::vector< T> < a class = "code hl_variable" href = "structgridfire_1_1_step_derivatives.html#ae0de268b86c2404379409c4feae0b34d" > dydt< / a > ; < / div >
< div class = "line" > < a id = "l00055" name = "l00055" > < / a > < span class = "lineno" > < a class = "line" href = "structgridfire_1_1_step_derivatives.html#ab4aeb41be952c7b5844e1ee81fef9008" > 55< / a > < / span > T < a class = "code hl_variable" href = "structgridfire_1_1_step_derivatives.html#ab4aeb41be952c7b5844e1ee81fef9008" > nuclearEnergyGenerationRate< / a > = T(0.0); < / div >
< div class = "line" > < a id = "l00056" name = "l00056" > < / a > < span class = "lineno" > 56< / span > };< / div >
2025-07-01 07:24:18 -04:00
< / div >
2025-07-01 11:40:51 -04:00
< div class = "line" > < a id = "l00057" name = "l00057" > < / a > < span class = "lineno" > 57< / span > < / div >
< div class = "foldopen" id = "foldopen00075" data-start = "{" data-end = "};" >
< div class = "line" > < a id = "l00075" name = "l00075" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_engine.html" > 75< / a > < / span > < span class = "keyword" > class < / span > < a class = "code hl_class" href = "classgridfire_1_1_engine.html" > Engine< / a > {< / div >
< div class = "line" > < a id = "l00076" name = "l00076" > < / a > < span class = "lineno" > 76< / span > < span class = "keyword" > public< / span > :< / div >
< div class = "line" > < a id = "l00080" name = "l00080" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_engine.html#a2e7970bed2100699f226f4141d5db037" > 80< / a > < / span > < span class = "keyword" > virtual< / span > < a class = "code hl_function" href = "classgridfire_1_1_engine.html#a2e7970bed2100699f226f4141d5db037" > ~Engine< / a > () = < span class = "keywordflow" > default< / span > ;< / div >
< div class = "line" > < a id = "l00081" name = "l00081" > < / a > < span class = "lineno" > 81< / span > < / div >
< div class = "line" > < a id = "l00086" name = "l00086" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_engine.html#a020e1b493d6964cafdad08fde697ceb3" > 86< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < span class = "keyword" > const< / span > std::vector< fourdst::atomic::Species> & < a class = "code hl_function" href = "classgridfire_1_1_engine.html#a020e1b493d6964cafdad08fde697ceb3" > getNetworkSpecies< / a > () < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00087" name = "l00087" > < / a > < span class = "lineno" > 87< / span > < / div >
< div class = "line" > < a id = "l00100" name = "l00100" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_engine.html#ac8b56124b6b49cd1802addb74a9a47c2" > 100< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < a class = "code hl_struct" href = "structgridfire_1_1_step_derivatives.html" > StepDerivatives< double> < / a > < a class = "code hl_function" href = "classgridfire_1_1_engine.html#ac8b56124b6b49cd1802addb74a9a47c2" > calculateRHSAndEnergy< / a > (< / div >
< div class = "line" > < a id = "l00101" name = "l00101" > < / a > < span class = "lineno" > 101< / span > < span class = "keyword" > const< / span > std::vector< double> & Y,< / div >
< div class = "line" > < a id = "l00102" name = "l00102" > < / a > < span class = "lineno" > 102< / span > < span class = "keywordtype" > double< / span > T9,< / div >
< div class = "line" > < a id = "l00103" name = "l00103" > < / a > < span class = "lineno" > 103< / span > < span class = "keywordtype" > double< / span > rho< / div >
< div class = "line" > < a id = "l00104" name = "l00104" > < / a > < span class = "lineno" > 104< / span > ) < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00105" name = "l00105" > < / a > < span class = "lineno" > 105< / span > };< / div >
2025-07-01 07:24:18 -04:00
< / div >
2025-07-01 11:40:51 -04:00
< div class = "line" > < a id = "l00106" name = "l00106" > < / a > < span class = "lineno" > 106< / span > < / div >
< div class = "foldopen" id = "foldopen00121" data-start = "{" data-end = "};" >
< div class = "line" > < a id = "l00121" name = "l00121" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html" > 121< / a > < / span > < span class = "keyword" > class < / span > < a class = "code hl_class" href = "classgridfire_1_1_dynamic_engine.html" > DynamicEngine< / a > : < span class = "keyword" > public< / span > < a class = "code hl_class" href = "classgridfire_1_1_engine.html" > Engine< / a > {< / div >
< div class = "line" > < a id = "l00122" name = "l00122" > < / a > < span class = "lineno" > 122< / span > < span class = "keyword" > public< / span > :< / div >
< div class = "line" > < a id = "l00133" name = "l00133" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a0a2fb3435ee3271ab9c806f225c61a7f" > 133< / a > < / span > < span class = "keyword" > virtual< / span > < span class = "keywordtype" > void< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a0a2fb3435ee3271ab9c806f225c61a7f" > generateJacobianMatrix< / a > (< / div >
< div class = "line" > < a id = "l00134" name = "l00134" > < / a > < span class = "lineno" > 134< / span > < span class = "keyword" > const< / span > std::vector< double> & Y,< / div >
< div class = "line" > < a id = "l00135" name = "l00135" > < / a > < span class = "lineno" > 135< / span > < span class = "keywordtype" > double< / span > T9, < span class = "keywordtype" > double< / span > rho< / div >
< div class = "line" > < a id = "l00136" name = "l00136" > < / a > < span class = "lineno" > 136< / span > ) = 0;< / div >
< div class = "line" > < a id = "l00137" name = "l00137" > < / a > < span class = "lineno" > 137< / span > < / div >
< div class = "line" > < a id = "l00147" name = "l00147" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a05d15ff35a6bc06a2fa7eda19838bd07" > 147< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < span class = "keywordtype" > double< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a05d15ff35a6bc06a2fa7eda19838bd07" > getJacobianMatrixEntry< / a > (< / div >
< div class = "line" > < a id = "l00148" name = "l00148" > < / a > < span class = "lineno" > 148< / span > < span class = "keywordtype" > int< / span > i,< / div >
< div class = "line" > < a id = "l00149" name = "l00149" > < / a > < span class = "lineno" > 149< / span > < span class = "keywordtype" > int< / span > j< / div >
< div class = "line" > < a id = "l00150" name = "l00150" > < / a > < span class = "lineno" > 150< / span > ) < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00151" name = "l00151" > < / a > < span class = "lineno" > 151< / span > < / div >
< div class = "line" > < a id = "l00158" name = "l00158" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#aeae6d84ef74d88fd2cdf07b82e98a16f" > 158< / a > < / span > < span class = "keyword" > virtual< / span > < span class = "keywordtype" > void< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#aeae6d84ef74d88fd2cdf07b82e98a16f" > generateStoichiometryMatrix< / a > () = 0;< / div >
< div class = "line" > < a id = "l00159" name = "l00159" > < / a > < span class = "lineno" > 159< / span > < / div >
< div class = "line" > < a id = "l00169" name = "l00169" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#afa108dd5227dbb1045e90d7b3bd8b84f" > 169< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < span class = "keywordtype" > int< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#afa108dd5227dbb1045e90d7b3bd8b84f" > getStoichiometryMatrixEntry< / a > (< / div >
< div class = "line" > < a id = "l00170" name = "l00170" > < / a > < span class = "lineno" > 170< / span > < span class = "keywordtype" > int< / span > speciesIndex,< / div >
< div class = "line" > < a id = "l00171" name = "l00171" > < / a > < span class = "lineno" > 171< / span > < span class = "keywordtype" > int< / span > reactionIndex< / div >
< div class = "line" > < a id = "l00172" name = "l00172" > < / a > < span class = "lineno" > 172< / span > ) < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00173" name = "l00173" > < / a > < span class = "lineno" > 173< / span > < / div >
< div class = "line" > < a id = "l00186" name = "l00186" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a6633b1757c41dd9e1c397333f4f9e785" > 186< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < span class = "keywordtype" > double< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a6633b1757c41dd9e1c397333f4f9e785" > calculateMolarReactionFlow< / a > (< / div >
< div class = "line" > < a id = "l00187" name = "l00187" > < / a > < span class = "lineno" > 187< / span > < span class = "keyword" > const< / span > < a class = "code hl_class" href = "classgridfire_1_1reaction_1_1_reaction.html" > reaction::Reaction< / a > & < a class = "code hl_namespace" href = "namespacegridfire_1_1reaction.html" > reaction< / a > ,< / div >
< div class = "line" > < a id = "l00188" name = "l00188" > < / a > < span class = "lineno" > 188< / span > < span class = "keyword" > const< / span > std::vector< double> & Y,< / div >
< div class = "line" > < a id = "l00189" name = "l00189" > < / a > < span class = "lineno" > 189< / span > < span class = "keywordtype" > double< / span > T9,< / div >
< div class = "line" > < a id = "l00190" name = "l00190" > < / a > < span class = "lineno" > 190< / span > < span class = "keywordtype" > double< / span > rho< / div >
< div class = "line" > < a id = "l00191" name = "l00191" > < / a > < span class = "lineno" > 191< / span > ) < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00192" name = "l00192" > < / a > < span class = "lineno" > 192< / span > < / div >
< div class = "line" > < a id = "l00198" name = "l00198" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#ad2a82099edbb374bbb2c9509ccdb1037" > 198< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < span class = "keyword" > const< / span > < a class = "code hl_typedef" href = "namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31" > reaction::LogicalReactionSet< / a > & < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#ad2a82099edbb374bbb2c9509ccdb1037" > getNetworkReactions< / a > () < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00199" name = "l00199" > < / a > < span class = "lineno" > 199< / span > < / div >
< div class = "line" > < a id = "l00211" name = "l00211" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a5d8ba98b230d2849035ee2507728fa15" > 211< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > std::unordered_map< fourdst::atomic::Species, double> < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a5d8ba98b230d2849035ee2507728fa15" > getSpeciesTimescales< / a > (< / div >
< div class = "line" > < a id = "l00212" name = "l00212" > < / a > < span class = "lineno" > 212< / span > < span class = "keyword" > const< / span > std::vector< double> & Y,< / div >
< div class = "line" > < a id = "l00213" name = "l00213" > < / a > < span class = "lineno" > 213< / span > < span class = "keywordtype" > double< / span > T9,< / div >
< div class = "line" > < a id = "l00214" name = "l00214" > < / a > < span class = "lineno" > 214< / span > < span class = "keywordtype" > double< / span > rho< / div >
< div class = "line" > < a id = "l00215" name = "l00215" > < / a > < span class = "lineno" > 215< / span > ) < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00216" name = "l00216" > < / a > < span class = "lineno" > 216< / span > < / div >
< div class = "line" > < a id = "l00217" name = "l00217" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#acd500e1cd788df1dc105d28a20dc5f4f" > 217< / a > < / span > < span class = "keyword" > virtual< / span > < span class = "keywordtype" > void< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#acd500e1cd788df1dc105d28a20dc5f4f" > update< / a > (< span class = "keyword" > const< / span > < a class = "code hl_struct" href = "structgridfire_1_1_net_in.html" > NetIn< / a > & netIn) = 0;< / div >
< div class = "line" > < a id = "l00218" name = "l00218" > < / a > < span class = "lineno" > 218< / span > < / div >
< div class = "line" > < a id = "l00219" name = "l00219" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a3fb44b6f55563a2f590f31916528f2bd" > 219< / a > < / span > < span class = "keyword" > virtual< / span > < span class = "keywordtype" > void< / span > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a3fb44b6f55563a2f590f31916528f2bd" > setScreeningModel< / a > (< a class = "code hl_enumeration" href = "namespacegridfire_1_1screening.html#aa82aafbc4f8c28d0a75b60798e3a7d25" > screening::ScreeningType< / a > model) = 0;< / div >
< div class = "line" > < a id = "l00220" name = "l00220" > < / a > < span class = "lineno" > 220< / span > < / div >
< div class = "line" > < a id = "l00221" name = "l00221" > < / a > < span class = "lineno" > < a class = "line" href = "classgridfire_1_1_dynamic_engine.html#a7a203f8e0f3a6744ddc912dfbcfdbcc0" > 221< / a > < / span > [[nodiscard]] < span class = "keyword" > virtual< / span > < a class = "code hl_enumeration" href = "namespacegridfire_1_1screening.html#aa82aafbc4f8c28d0a75b60798e3a7d25" > screening::ScreeningType< / a > < a class = "code hl_function" href = "classgridfire_1_1_dynamic_engine.html#a7a203f8e0f3a6744ddc912dfbcfdbcc0" > getScreeningModel< / a > () < span class = "keyword" > const< / span > = 0;< / div >
< div class = "line" > < a id = "l00222" name = "l00222" > < / a > < span class = "lineno" > 222< / span > };< / div >
2025-07-01 07:24:18 -04:00
< / div >
2025-07-01 11:40:51 -04:00
< div class = "line" > < a id = "l00223" name = "l00223" > < / a > < span class = "lineno" > 223< / span > }< / div >
2025-07-01 07:24:18 -04:00
< / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html" > gridfire::DynamicEngine< / a > < / div > < div class = "ttdoc" > Abstract class for engines supporting Jacobian and stoichiometry operations.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00121" > engine_abstract.h:121< / a > < / div > < / div >
2025-07-01 07:24:18 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a05d15ff35a6bc06a2fa7eda19838bd07" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a05d15ff35a6bc06a2fa7eda19838bd07" > gridfire::DynamicEngine::getJacobianMatrixEntry< / a > < / div > < div class = "ttdeci" > virtual double getJacobianMatrixEntry(int i, int j) const =0< / div > < div class = "ttdoc" > Get an entry from the previously generated Jacobian matrix.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a0a2fb3435ee3271ab9c806f225c61a7f" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a0a2fb3435ee3271ab9c806f225c61a7f" > gridfire::DynamicEngine::generateJacobianMatrix< / a > < / div > < div class = "ttdeci" > virtual void generateJacobianMatrix(const std::vector< double > & Y, double T9, double rho)=0< / div > < div class = "ttdoc" > Generate the Jacobian matrix for the current state.< / div > < / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a3fb44b6f55563a2f590f31916528f2bd" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a3fb44b6f55563a2f590f31916528f2bd" > gridfire::DynamicEngine::setScreeningModel< / a > < / div > < div class = "ttdeci" > virtual void setScreeningModel(screening::ScreeningType model)=0< / div > < / div >
2025-07-01 07:24:18 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a5d8ba98b230d2849035ee2507728fa15" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a5d8ba98b230d2849035ee2507728fa15" > gridfire::DynamicEngine::getSpeciesTimescales< / a > < / div > < div class = "ttdeci" > virtual std::unordered_map< fourdst::atomic::Species, double > getSpeciesTimescales(const std::vector< double > & Y, double T9, double rho) const =0< / div > < div class = "ttdoc" > Compute timescales for all species in the network.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a6633b1757c41dd9e1c397333f4f9e785" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a6633b1757c41dd9e1c397333f4f9e785" > gridfire::DynamicEngine::calculateMolarReactionFlow< / a > < / div > < div class = "ttdeci" > virtual double calculateMolarReactionFlow(const reaction::Reaction & reaction, const std::vector< double > & Y, double T9, double rho) const =0< / div > < div class = "ttdoc" > Calculate the molar reaction flow for a given reaction.< / div > < / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_a7a203f8e0f3a6744ddc912dfbcfdbcc0" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#a7a203f8e0f3a6744ddc912dfbcfdbcc0" > gridfire::DynamicEngine::getScreeningModel< / a > < / div > < div class = "ttdeci" > virtual screening::ScreeningType getScreeningModel() const =0< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_acd500e1cd788df1dc105d28a20dc5f4f" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#acd500e1cd788df1dc105d28a20dc5f4f" > gridfire::DynamicEngine::update< / a > < / div > < div class = "ttdeci" > virtual void update(const NetIn & netIn)=0< / div > < / div >
2025-07-01 07:24:18 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_ad2a82099edbb374bbb2c9509ccdb1037" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#ad2a82099edbb374bbb2c9509ccdb1037" > gridfire::DynamicEngine::getNetworkReactions< / a > < / div > < div class = "ttdeci" > virtual const reaction::LogicalReactionSet & getNetworkReactions() const =0< / div > < div class = "ttdoc" > Get the set of logical reactions in the network.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_aeae6d84ef74d88fd2cdf07b82e98a16f" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#aeae6d84ef74d88fd2cdf07b82e98a16f" > gridfire::DynamicEngine::generateStoichiometryMatrix< / a > < / div > < div class = "ttdeci" > virtual void generateStoichiometryMatrix()=0< / div > < div class = "ttdoc" > Generate the stoichiometry matrix for the network.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_dynamic_engine_html_afa108dd5227dbb1045e90d7b3bd8b84f" > < div class = "ttname" > < a href = "classgridfire_1_1_dynamic_engine.html#afa108dd5227dbb1045e90d7b3bd8b84f" > gridfire::DynamicEngine::getStoichiometryMatrixEntry< / a > < / div > < div class = "ttdeci" > virtual int getStoichiometryMatrixEntry(int speciesIndex, int reactionIndex) const =0< / div > < div class = "ttdoc" > Get an entry from the stoichiometry matrix.< / div > < / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_engine_html" > < div class = "ttname" > < a href = "classgridfire_1_1_engine.html" > gridfire::Engine< / a > < / div > < div class = "ttdoc" > Abstract base class for a reaction network engine.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00075" > engine_abstract.h:75< / a > < / div > < / div >
2025-07-01 07:24:18 -04:00
< div class = "ttc" id = "aclassgridfire_1_1_engine_html_a020e1b493d6964cafdad08fde697ceb3" > < div class = "ttname" > < a href = "classgridfire_1_1_engine.html#a020e1b493d6964cafdad08fde697ceb3" > gridfire::Engine::getNetworkSpecies< / a > < / div > < div class = "ttdeci" > virtual const std::vector< fourdst::atomic::Species > & getNetworkSpecies() const =0< / div > < div class = "ttdoc" > Get the list of species in the network.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_engine_html_a2e7970bed2100699f226f4141d5db037" > < div class = "ttname" > < a href = "classgridfire_1_1_engine.html#a2e7970bed2100699f226f4141d5db037" > gridfire::Engine::~Engine< / a > < / div > < div class = "ttdeci" > virtual ~Engine()=default< / div > < div class = "ttdoc" > Virtual destructor.< / div > < / div >
< div class = "ttc" id = "aclassgridfire_1_1_engine_html_ac8b56124b6b49cd1802addb74a9a47c2" > < div class = "ttname" > < a href = "classgridfire_1_1_engine.html#ac8b56124b6b49cd1802addb74a9a47c2" > gridfire::Engine::calculateRHSAndEnergy< / a > < / div > < div class = "ttdeci" > virtual StepDerivatives< double > calculateRHSAndEnergy(const std::vector< double > & Y, double T9, double rho) const =0< / div > < div class = "ttdoc" > Calculate the right-hand side (dY/dt) and energy generation.< / div > < / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "aclassgridfire_1_1reaction_1_1_reaction_html" > < div class = "ttname" > < a href = "classgridfire_1_1reaction_1_1_reaction.html" > gridfire::reaction::Reaction< / a > < / div > < div class = "ttdoc" > Represents a single nuclear reaction from a specific data source.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "reaction_8h_source.html#l00072" > reaction.h:72< / a > < / div > < / div >
< div class = "ttc" id = "aconceptgridfire_1_1_is_arithmetic_or_a_d_html" > < div class = "ttname" > < a href = "conceptgridfire_1_1_is_arithmetic_or_a_d.html" > gridfire::IsArithmeticOrAD< / a > < / div > < div class = "ttdoc" > Concept for types allowed in engine calculations.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00033" > engine_abstract.h:33< / a > < / div > < / div >
< div class = "ttc" id = "anamespacegridfire_1_1reaction_html" > < div class = "ttname" > < a href = "namespacegridfire_1_1reaction.html" > gridfire::reaction< / a > < / div > < div class = "ttdef" > < b > Definition< / b > < a href = "reaction_8h_source.html#l00025" > reaction.h:25< / a > < / div > < / div >
< div class = "ttc" id = "anamespacegridfire_1_1reaction_html_aa86f08712565f278adacc7cd2361eb31" > < div class = "ttname" > < a href = "namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31" > gridfire::reaction::LogicalReactionSet< / a > < / div > < div class = "ttdeci" > TemplatedReactionSet< LogicalReaction > LogicalReactionSet< / div > < div class = "ttdoc" > A set of logical reactions.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "reaction_8h_source.html#l00557" > reaction.h:557< / a > < / div > < / div >
< div class = "ttc" id = "anamespacegridfire_1_1screening_html_aa82aafbc4f8c28d0a75b60798e3a7d25" > < div class = "ttname" > < a href = "namespacegridfire_1_1screening.html#aa82aafbc4f8c28d0a75b60798e3a7d25" > gridfire::screening::ScreeningType< / a > < / div > < div class = "ttdeci" > ScreeningType< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "screening__types_8h_source.html#l00008" > screening_types.h:8< / a > < / div > < / div >
< div class = "ttc" id = "anamespacegridfire_html" > < div class = "ttname" > < a href = "namespacegridfire.html" > gridfire< / a > < / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00024" > engine_abstract.h:24< / a > < / div > < / div >
< div class = "ttc" id = "anetwork_8h_html" > < div class = "ttname" > < a href = "network_8h.html" > network.h< / a > < / div > < / div >
2025-07-01 07:24:18 -04:00
< div class = "ttc" id = "areaction_8h_html" > < div class = "ttname" > < a href = "reaction_8h.html" > reaction.h< / a > < / div > < div class = "ttdoc" > Defines classes for representing and managing nuclear reactions.< / div > < / div >
2025-07-01 11:40:51 -04:00
< div class = "ttc" id = "ascreening__abstract_8h_html" > < div class = "ttname" > < a href = "screening__abstract_8h.html" > screening_abstract.h< / a > < / div > < / div >
< div class = "ttc" id = "ascreening__types_8h_html" > < div class = "ttname" > < a href = "screening__types_8h.html" > screening_types.h< / a > < / div > < / div >
< div class = "ttc" id = "astructgridfire_1_1_net_in_html" > < div class = "ttname" > < a href = "structgridfire_1_1_net_in.html" > gridfire::NetIn< / a > < / div > < div class = "ttdef" > < b > Definition< / b > < a href = "network_8h_source.html#l00053" > network.h:53< / a > < / div > < / div >
< div class = "ttc" id = "astructgridfire_1_1_step_derivatives_html" > < div class = "ttname" > < a href = "structgridfire_1_1_step_derivatives.html" > gridfire::StepDerivatives< / a > < / div > < div class = "ttdoc" > Structure holding derivatives and energy generation for a network step.< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00053" > engine_abstract.h:53< / a > < / div > < / div >
< div class = "ttc" id = "astructgridfire_1_1_step_derivatives_html_ab4aeb41be952c7b5844e1ee81fef9008" > < div class = "ttname" > < a href = "structgridfire_1_1_step_derivatives.html#ab4aeb41be952c7b5844e1ee81fef9008" > gridfire::StepDerivatives::nuclearEnergyGenerationRate< / a > < / div > < div class = "ttdeci" > T nuclearEnergyGenerationRate< / div > < div class = "ttdoc" > Specific energy generation rate (e.g., erg/g/s).< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00055" > engine_abstract.h:55< / a > < / div > < / div >
< div class = "ttc" id = "astructgridfire_1_1_step_derivatives_html_ae0de268b86c2404379409c4feae0b34d" > < div class = "ttname" > < a href = "structgridfire_1_1_step_derivatives.html#ae0de268b86c2404379409c4feae0b34d" > gridfire::StepDerivatives::dydt< / a > < / div > < div class = "ttdeci" > std::vector< T > dydt< / div > < div class = "ttdoc" > Derivatives of abundances (dY/dt for each species).< / div > < div class = "ttdef" > < b > Definition< / b > < a href = "#l00054" > engine_abstract.h:54< / a > < / div > < / div >
2025-07-01 07:24:18 -04:00
< / div > <!-- fragment --> < / div > <!-- contents -->
< / div > <!-- doc - content -->
<!-- start footer part -->
< div id = "nav-path" class = "navpath" > <!-- id is needed for treeview function! -->
< ul >
< li class = "navelem" > < a class = "el" href = "dir_68267d1309a1af8e8297ef4c3efbcdba.html" > src< / a > < / li > < li class = "navelem" > < a class = "el" href = "dir_fc4c7f03e1a69a98c370fae55a743828.html" > network< / a > < / li > < li class = "navelem" > < a class = "el" href = "dir_5cccfa813acdf3744b542715860d37b2.html" > include< / a > < / li > < li class = "navelem" > < a class = "el" href = "dir_2a1262ef5950eb718393488a3eb5aa9f.html" > gridfire< / a > < / li > < li class = "navelem" > < a class = "el" href = "dir_6b2e1e22dfdea3280d50981209bf7529.html" > engine< / a > < / li > < li class = "navelem" > < a class = "el" href = "engine__abstract_8h.html" > engine_abstract.h< / a > < / li >
< li class = "footer" > Generated by < a href = "https://www.doxygen.org/index.html" > < img class = "footer" src = "doxygen.svg" width = "104" height = "31" alt = "doxygen" / > < / a > 1.13.2 < / li >
< / ul >
< / div >
< / body >
< / html >