Files
GridFire/docs/html/engine__adaptive_8h_source.html

293 lines
56 KiB
HTML

<!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/views/engine_adaptive.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">&#160;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&amp;dn=expat.txt MIT */
var searchBox = new SearchBox("searchBox", "search/",'.html');
/* @license-end */
</script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&amp;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&amp;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&amp;dn=expat.txt MIT */
$(function(){initNavTree('engine__adaptive_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_adaptive.h</div></div>
</div><!--header-->
<div class="contents">
<a href="engine__adaptive_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><span class="preprocessor">#include &quot;<a class="code" href="engine__abstract_8h.html">gridfire/engine/engine_abstract.h</a>&quot;</span></div>
<div class="line"><a id="l00003" name="l00003"></a><span class="lineno"> 3</span><span class="preprocessor">#include &quot;<a class="code" href="engine__view__abstract_8h.html">gridfire/engine/views/engine_view_abstract.h</a>&quot;</span></div>
<div class="line"><a id="l00004" name="l00004"></a><span class="lineno"> 4</span><span class="preprocessor">#include &quot;<a class="code" href="screening__abstract_8h.html">gridfire/screening/screening_abstract.h</a>&quot;</span></div>
<div class="line"><a id="l00005" name="l00005"></a><span class="lineno"> 5</span><span class="preprocessor">#include &quot;<a class="code" href="screening__types_8h.html">gridfire/screening/screening_types.h</a>&quot;</span></div>
<div class="line"><a id="l00006" name="l00006"></a><span class="lineno"> 6</span><span class="preprocessor">#include &quot;<a class="code" href="network_8h.html">gridfire/network.h</a>&quot;</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 &quot;fourdst/composition/atomicSpecies.h&quot;</span></div>
<div class="line"><a id="l00009" name="l00009"></a><span class="lineno"> 9</span><span class="preprocessor">#include &quot;fourdst/config/config.h&quot;</span></div>
<div class="line"><a id="l00010" name="l00010"></a><span class="lineno"> 10</span><span class="preprocessor">#include &quot;fourdst/logging/logging.h&quot;</span></div>
<div class="line"><a id="l00011" name="l00011"></a><span class="lineno"> 11</span> </div>
<div class="line"><a id="l00012" name="l00012"></a><span class="lineno"> 12</span><span class="preprocessor">#include &quot;quill/Logger.h&quot;</span></div>
<div class="line"><a id="l00013" name="l00013"></a><span class="lineno"> 13</span> </div>
<div class="line"><a id="l00014" name="l00014"></a><span class="lineno"> 14</span><span class="keyword">namespace </span><a class="code hl_namespace" href="namespacegridfire.html">gridfire</a> {</div>
<div class="foldopen" id="foldopen00047" data-start="{" data-end="};">
<div class="line"><a id="l00047" name="l00047"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html"> 47</a></span> <span class="keyword">class </span><a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#ad599363cdd457e72e2e2937b0222c455">AdaptiveEngineView</a> final : <span class="keyword">public</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_view.html">EngineView</a>&lt;DynamicEngine&gt; {</div>
<div class="line"><a id="l00048" name="l00048"></a><span class="lineno"> 48</span> <span class="keyword">public</span>:</div>
<div class="line"><a id="l00057" name="l00057"></a><span class="lineno"> 57</span> <span class="keyword">explicit</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#ad599363cdd457e72e2e2937b0222c455">AdaptiveEngineView</a>(<a class="code hl_class" href="classgridfire_1_1_dynamic_engine.html">DynamicEngine</a>&amp; baseEngine);</div>
<div class="line"><a id="l00058" name="l00058"></a><span class="lineno"> 58</span></div>
<div class="line"><a id="l00076" name="l00076"></a><span class="lineno"> 76</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a7d61e73f5158f1574cda3edc90c51f7e">update</a>(<span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>&amp; netIn) <span class="keyword">override</span>;</div>
<div class="line"><a id="l00077" name="l00077"></a><span class="lineno"> 77</span></div>
<div class="line"><a id="l00082" name="l00082"></a><span class="lineno"> 82</span> [[nodiscard]] <span class="keyword">const</span> std::vector&lt;fourdst::atomic::Species&gt;&amp; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#ac83a8efe25c0e5b9bf7756ac3a500bb1">getNetworkSpecies</a>() <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00083" name="l00083"></a><span class="lineno"> 83</span></div>
<div class="line"><a id="l00100" name="l00100"></a><span class="lineno"> 100</span> [[nodiscard]] <a class="code hl_struct" href="structgridfire_1_1_step_derivatives.html">StepDerivatives&lt;double&gt;</a> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a7b276b7210be588263395bdb0497fc6d">calculateRHSAndEnergy</a>(</div>
<div class="line"><a id="l00101" name="l00101"></a><span class="lineno"> 101</span> <span class="keyword">const</span> std::vector&lt;double&gt; &amp;Y_culled,</div>
<div class="line"><a id="l00102" name="l00102"></a><span class="lineno"> 102</span> <span class="keyword">const</span> <span class="keywordtype">double</span> T9,</div>
<div class="line"><a id="l00103" name="l00103"></a><span class="lineno"> 103</span> <span class="keyword">const</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 override</span>;</div>
<div class="line"><a id="l00105" name="l00105"></a><span class="lineno"> 105</span></div>
<div class="line"><a id="l00119" name="l00119"></a><span class="lineno"> 119</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#ac9aab6f60e80a9228b2b19b1b10449ef">generateJacobianMatrix</a>(</div>
<div class="line"><a id="l00120" name="l00120"></a><span class="lineno"> 120</span> <span class="keyword">const</span> std::vector&lt;double&gt; &amp;Y_culled,</div>
<div class="line"><a id="l00121" name="l00121"></a><span class="lineno"> 121</span> <span class="keyword">const</span> <span class="keywordtype">double</span> T9,</div>
<div class="line"><a id="l00122" name="l00122"></a><span class="lineno"> 122</span> <span class="keyword">const</span> <span class="keywordtype">double</span> rho</div>
<div class="line"><a id="l00123" name="l00123"></a><span class="lineno"> 123</span> ) <span class="keyword">override</span>;</div>
<div class="line"><a id="l00124" name="l00124"></a><span class="lineno"> 124</span></div>
<div class="line"><a id="l00139" name="l00139"></a><span class="lineno"> 139</span> [[nodiscard]] <span class="keywordtype">double</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a4710d218c8a0fd161e994ecd60b48e58">getJacobianMatrixEntry</a>(</div>
<div class="line"><a id="l00140" name="l00140"></a><span class="lineno"> 140</span> <span class="keyword">const</span> <span class="keywordtype">int</span> i_culled,</div>
<div class="line"><a id="l00141" name="l00141"></a><span class="lineno"> 141</span> <span class="keyword">const</span> <span class="keywordtype">int</span> j_culled</div>
<div class="line"><a id="l00142" name="l00142"></a><span class="lineno"> 142</span> ) <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00143" name="l00143"></a><span class="lineno"> 143</span></div>
<div class="line"><a id="l00153" name="l00153"></a><span class="lineno"> 153</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a231193a61ba5a31e8eb92b0d4ce69111">generateStoichiometryMatrix</a>() <span class="keyword">override</span>;</div>
<div class="line"><a id="l00154" name="l00154"></a><span class="lineno"> 154</span></div>
<div class="line"><a id="l00169" name="l00169"></a><span class="lineno"> 169</span> [[nodiscard]] <span class="keywordtype">int</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a67b4ea8cad115394bb4a42cc39d696f9">getStoichiometryMatrixEntry</a>(</div>
<div class="line"><a id="l00170" name="l00170"></a><span class="lineno"> 170</span> <span class="keyword">const</span> <span class="keywordtype">int</span> speciesIndex_culled,</div>
<div class="line"><a id="l00171" name="l00171"></a><span class="lineno"> 171</span> <span class="keyword">const</span> <span class="keywordtype">int</span> reactionIndex_culled</div>
<div class="line"><a id="l00172" name="l00172"></a><span class="lineno"> 172</span> ) <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00173" name="l00173"></a><span class="lineno"> 173</span></div>
<div class="line"><a id="l00189" name="l00189"></a><span class="lineno"> 189</span> [[nodiscard]] <span class="keywordtype">double</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a048d4b1d41ecb4125a558d1b9ed7cb31">calculateMolarReactionFlow</a>(</div>
<div class="line"><a id="l00190" name="l00190"></a><span class="lineno"> 190</span> <span class="keyword">const</span> <a class="code hl_class" href="classgridfire_1_1reaction_1_1_reaction.html">reaction::Reaction</a> &amp;<a class="code hl_namespace" href="namespacegridfire_1_1reaction.html">reaction</a>,</div>
<div class="line"><a id="l00191" name="l00191"></a><span class="lineno"> 191</span> <span class="keyword">const</span> std::vector&lt;double&gt; &amp;Y_culled,</div>
<div class="line"><a id="l00192" name="l00192"></a><span class="lineno"> 192</span> <span class="keywordtype">double</span> T9,</div>
<div class="line"><a id="l00193" name="l00193"></a><span class="lineno"> 193</span> <span class="keywordtype">double</span> rho</div>
<div class="line"><a id="l00194" name="l00194"></a><span class="lineno"> 194</span> ) <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00195" name="l00195"></a><span class="lineno"> 195</span></div>
<div class="line"><a id="l00201" name="l00201"></a><span class="lineno"> 201</span> [[nodiscard]] <span class="keyword">const</span> <a class="code hl_typedef" href="namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31">reaction::LogicalReactionSet</a>&amp; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a12cc2f352678fba9688363ba1876ab9c">getNetworkReactions</a>() <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00202" name="l00202"></a><span class="lineno"> 202</span></div>
<div class="line"><a id="l00216" name="l00216"></a><span class="lineno"> 216</span> [[nodiscard]] std::unordered_map&lt;fourdst::atomic::Species, double&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a93b38d0fdc4647f6f7340172dae17872">getSpeciesTimescales</a>(</div>
<div class="line"><a id="l00217" name="l00217"></a><span class="lineno"> 217</span> <span class="keyword">const</span> std::vector&lt;double&gt; &amp;Y_culled,</div>
<div class="line"><a id="l00218" name="l00218"></a><span class="lineno"> 218</span> <span class="keywordtype">double</span> T9,</div>
<div class="line"><a id="l00219" name="l00219"></a><span class="lineno"> 219</span> <span class="keywordtype">double</span> rho</div>
<div class="line"><a id="l00220" name="l00220"></a><span class="lineno"> 220</span> ) <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00221" name="l00221"></a><span class="lineno"> 221</span></div>
<div class="line"><a id="l00226" name="l00226"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#aee095b30a9dce5fcb5ae2fa1d2aa192c"> 226</a></span> [[nodiscard]] <span class="keyword">const</span> <a class="code hl_class" href="classgridfire_1_1_dynamic_engine.html">DynamicEngine</a>&amp; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#aee095b30a9dce5fcb5ae2fa1d2aa192c">getBaseEngine</a>()<span class="keyword"> const override </span>{ <span class="keywordflow">return</span> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8">m_baseEngine</a>; }</div>
<div class="line"><a id="l00227" name="l00227"></a><span class="lineno"> 227</span></div>
<div class="line"><a id="l00243" name="l00243"></a><span class="lineno"> 243</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#aae4ddbef1c4e2202fd236221a4bf376b">setScreeningModel</a>(<a class="code hl_enumeration" href="namespacegridfire_1_1screening.html#aa82aafbc4f8c28d0a75b60798e3a7d25">screening::ScreeningType</a> model) <span class="keyword">override</span>;</div>
<div class="line"><a id="l00244" name="l00244"></a><span class="lineno"> 244</span></div>
<div class="line"><a id="l00258" name="l00258"></a><span class="lineno"> 258</span> [[nodiscard]] <a class="code hl_enumeration" href="namespacegridfire_1_1screening.html#aa82aafbc4f8c28d0a75b60798e3a7d25">screening::ScreeningType</a> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a0ab1199f900a58f309c3c36532c9164f">getScreeningModel</a>() <span class="keyword">const override</span>;</div>
<div class="line"><a id="l00259" name="l00259"></a><span class="lineno"> 259</span> <span class="keyword">private</span>:</div>
<div class="line"><a id="l00260" name="l00260"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3"> 260</a></span> <span class="keyword">using </span><a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3">Config</a> = fourdst::config::Config;</div>
<div class="line"><a id="l00261" name="l00261"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a5eaf7c3a4e28cd3a4f34979b88a80103"> 261</a></span> <span class="keyword">using </span><a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#a5eaf7c3a4e28cd3a4f34979b88a80103">LogManager</a> = fourdst::logging::LogManager;</div>
<div class="line"><a id="l00263" name="l00263"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a14171a9ccc45a63996a967c72983de30"> 263</a></span> <a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3">Config</a>&amp; <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a14171a9ccc45a63996a967c72983de30">m_config</a> = Config::getInstance();</div>
<div class="line"><a id="l00265" name="l00265"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#ac5bdbe46f87d38d9f23ece5743dcd193"> 265</a></span> quill::Logger* <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#ac5bdbe46f87d38d9f23ece5743dcd193">m_logger</a> = LogManager::getInstance().getLogger(<span class="stringliteral">&quot;log&quot;</span>);</div>
<div class="line"><a id="l00266" name="l00266"></a><span class="lineno"> 266</span></div>
<div class="line"><a id="l00268" name="l00268"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8"> 268</a></span> <a class="code hl_class" href="classgridfire_1_1_dynamic_engine.html">DynamicEngine</a>&amp; <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8">m_baseEngine</a>;</div>
<div class="line"><a id="l00269" name="l00269"></a><span class="lineno"> 269</span></div>
<div class="line"><a id="l00271" name="l00271"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#af4bc9fc6e4afcd6a53c49ca6e2a95940"> 271</a></span> std::vector&lt;fourdst::atomic::Species&gt; <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#af4bc9fc6e4afcd6a53c49ca6e2a95940">m_activeSpecies</a>;</div>
<div class="line"><a id="l00273" name="l00273"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a19fc7e02e216b797aa643fa35e429800"> 273</a></span> <a class="code hl_typedef" href="namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31">reaction::LogicalReactionSet</a> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a19fc7e02e216b797aa643fa35e429800">m_activeReactions</a>;</div>
<div class="line"><a id="l00274" name="l00274"></a><span class="lineno"> 274</span></div>
<div class="line"><a id="l00276" name="l00276"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a5f66204a0ff5b27eed243afddecb0093"> 276</a></span> std::vector&lt;size_t&gt; <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a5f66204a0ff5b27eed243afddecb0093">m_speciesIndexMap</a>;</div>
<div class="line"><a id="l00278" name="l00278"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a21c6e33bbf8c18fd5b5eaabb469054de"> 278</a></span> std::vector&lt;size_t&gt; <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a21c6e33bbf8c18fd5b5eaabb469054de">m_reactionIndexMap</a>;</div>
<div class="line"><a id="l00279" name="l00279"></a><span class="lineno"> 279</span></div>
<div class="line"><a id="l00281" name="l00281"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a63580db57e0f48f508906a11ccfd465e"> 281</a></span> <span class="keywordtype">bool</span> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a63580db57e0f48f508906a11ccfd465e">m_isStale</a> = <span class="keyword">true</span>;</div>
<div class="line"><a id="l00282" name="l00282"></a><span class="lineno"> 282</span> </div>
<div class="line"><a id="l00283" name="l00283"></a><span class="lineno"> 283</span> <span class="keyword">private</span>:</div>
<div class="foldopen" id="foldopen00287" data-start="{" data-end="};">
<div class="line"><a id="l00287" name="l00287"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html"> 287</a></span> <span class="keyword">struct </span><a class="code hl_struct" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html">ReactionFlow</a> {</div>
<div class="line"><a id="l00288" name="l00288"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a3bb21f20df8115d37108cf3c3be3bc6f"> 288</a></span> <span class="keyword">const</span> <a class="code hl_class" href="classgridfire_1_1reaction_1_1_logical_reaction.html">reaction::LogicalReaction</a>* <a class="code hl_variable" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a3bb21f20df8115d37108cf3c3be3bc6f">reactionPtr</a>;</div>
<div class="line"><a id="l00289" name="l00289"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a6c2e9087f6c8af5d89a5f0de7bd4a5b4"> 289</a></span> <span class="keywordtype">double</span> <a class="code hl_variable" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a6c2e9087f6c8af5d89a5f0de7bd4a5b4">flowRate</a>;</div>
<div class="line"><a id="l00290" name="l00290"></a><span class="lineno"> 290</span> };</div>
</div>
<div class="line"><a id="l00291" name="l00291"></a><span class="lineno"> 291</span> <span class="keyword">private</span>:</div>
<div class="line"><a id="l00302" name="l00302"></a><span class="lineno"> 302</span> [[nodiscard]] std::vector&lt;size_t&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a896d29325b4233e83d9298850b617a2d">constructSpeciesIndexMap</a>() <span class="keyword">const</span>;</div>
<div class="line"><a id="l00303" name="l00303"></a><span class="lineno"> 303</span></div>
<div class="line"><a id="l00314" name="l00314"></a><span class="lineno"> 314</span> [[nodiscard]] std::vector&lt;size_t&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a89614f4a48f60c4170a0197f45303e7c">constructReactionIndexMap</a>() <span class="keyword">const</span>;</div>
<div class="line"><a id="l00315" name="l00315"></a><span class="lineno"> 315</span></div>
<div class="line"><a id="l00323" name="l00323"></a><span class="lineno"> 323</span> [[nodiscard]] std::vector&lt;double&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a68695f056b660e91285b7e5a931612e1">mapCulledToFull</a>(<span class="keyword">const</span> std::vector&lt;double&gt;&amp; culled) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00324" name="l00324"></a><span class="lineno"> 324</span></div>
<div class="line"><a id="l00332" name="l00332"></a><span class="lineno"> 332</span> [[nodiscard]] std::vector&lt;double&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a3d9d8e862d1c2f0a8ba460c57f6a7f44">mapFullToCulled</a>(<span class="keyword">const</span> std::vector&lt;double&gt;&amp; full) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00333" name="l00333"></a><span class="lineno"> 333</span></div>
<div class="line"><a id="l00342" name="l00342"></a><span class="lineno"> 342</span> [[nodiscard]] <span class="keywordtype">size_t</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a256d14a333f9401039b826cc889761a8">mapCulledToFullSpeciesIndex</a>(<span class="keywordtype">size_t</span> culledSpeciesIndex) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00343" name="l00343"></a><span class="lineno"> 343</span></div>
<div class="line"><a id="l00352" name="l00352"></a><span class="lineno"> 352</span> [[nodiscard]] <span class="keywordtype">size_t</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a91e742642d8a8d9ec0620779927e5101">mapCulledToFullReactionIndex</a>(<span class="keywordtype">size_t</span> culledReactionIndex) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00353" name="l00353"></a><span class="lineno"> 353</span></div>
<div class="line"><a id="l00359" name="l00359"></a><span class="lineno"> 359</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#aedc0dedb51c81e03f253cc409a5d5c40">validateState</a>() <span class="keyword">const</span>;</div>
<div class="line"><a id="l00360" name="l00360"></a><span class="lineno"> 360</span></div>
<div class="line"><a id="l00382" name="l00382"></a><span class="lineno"> 382</span> std::vector&lt;ReactionFlow&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#abdbaf4b87629efe43ac1255dad424c0c">calculateAllReactionFlows</a>(</div>
<div class="line"><a id="l00383" name="l00383"></a><span class="lineno"> 383</span> <span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>&amp; netIn,</div>
<div class="line"><a id="l00384" name="l00384"></a><span class="lineno"> 384</span> std::vector&lt;double&gt;&amp; out_Y_Full</div>
<div class="line"><a id="l00385" name="l00385"></a><span class="lineno"> 385</span> ) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00403" name="l00403"></a><span class="lineno"> 403</span> [[nodiscard]] std::unordered_set&lt;fourdst::atomic::Species&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a0ed21f7e7c1034fc87b40d4116c4221b">findReachableSpecies</a>(</div>
<div class="line"><a id="l00404" name="l00404"></a><span class="lineno"> 404</span> <span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>&amp; netIn</div>
<div class="line"><a id="l00405" name="l00405"></a><span class="lineno"> 405</span> ) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00426" name="l00426"></a><span class="lineno"> 426</span> [[nodiscard]] std::vector&lt;const reaction::LogicalReaction*&gt; <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a42417e96fe9fd623458af109401daf08">cullReactionsByFlow</a>(</div>
<div class="line"><a id="l00427" name="l00427"></a><span class="lineno"> 427</span> <span class="keyword">const</span> std::vector&lt;ReactionFlow&gt;&amp; allFlows,</div>
<div class="line"><a id="l00428" name="l00428"></a><span class="lineno"> 428</span> <span class="keyword">const</span> std::unordered_set&lt;fourdst::atomic::Species&gt;&amp; reachableSpecies,</div>
<div class="line"><a id="l00429" name="l00429"></a><span class="lineno"> 429</span> <span class="keyword">const</span> std::vector&lt;double&gt;&amp; Y_full,</div>
<div class="line"><a id="l00430" name="l00430"></a><span class="lineno"> 430</span> <span class="keywordtype">double</span> maxFlow</div>
<div class="line"><a id="l00431" name="l00431"></a><span class="lineno"> 431</span> ) <span class="keyword">const</span>;</div>
<div class="line"><a id="l00447" name="l00447"></a><span class="lineno"> 447</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#aa79fb382c98461b02a2c30668491e6c5">finalizeActiveSet</a>(</div>
<div class="line"><a id="l00448" name="l00448"></a><span class="lineno"> 448</span> <span class="keyword">const</span> std::vector&lt;const reaction::LogicalReaction*&gt;&amp; finalReactions</div>
<div class="line"><a id="l00449" name="l00449"></a><span class="lineno"> 449</span> );</div>
<div class="line"><a id="l00450" name="l00450"></a><span class="lineno"> 450</span> };</div>
</div>
<div class="line"><a id="l00451" name="l00451"></a><span class="lineno"> 451</span>}</div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a048d4b1d41ecb4125a558d1b9ed7cb31"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a048d4b1d41ecb4125a558d1b9ed7cb31">gridfire::AdaptiveEngineView::calculateMolarReactionFlow</a></div><div class="ttdeci">double calculateMolarReactionFlow(const reaction::Reaction &amp;reaction, const std::vector&lt; double &gt; &amp;Y_culled, double T9, double rho) const override</div><div class="ttdoc">Calculates the molar reaction flow for a given reaction in the active network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00175">engine_adaptive.cpp:175</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a0ab1199f900a58f309c3c36532c9164f"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a0ab1199f900a58f309c3c36532c9164f">gridfire::AdaptiveEngineView::getScreeningModel</a></div><div class="ttdeci">screening::ScreeningType getScreeningModel() const override</div><div class="ttdoc">Gets the screening model from the base engine.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00220">engine_adaptive.cpp:220</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a0ed21f7e7c1034fc87b40d4116c4221b"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a0ed21f7e7c1034fc87b40d4116c4221b">gridfire::AdaptiveEngineView::findReachableSpecies</a></div><div class="ttdeci">std::unordered_set&lt; fourdst::atomic::Species &gt; findReachableSpecies(const NetIn &amp;netIn) const</div><div class="ttdoc">Finds all species that are reachable from the initial fuel through the reaction network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00299">engine_adaptive.cpp:299</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a12cc2f352678fba9688363ba1876ab9c"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a12cc2f352678fba9688363ba1876ab9c">gridfire::AdaptiveEngineView::getNetworkReactions</a></div><div class="ttdeci">const reaction::LogicalReactionSet &amp; getNetworkReactions() const override</div><div class="ttdoc">Gets the set of active logical reactions in the network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00192">engine_adaptive.cpp:192</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a14171a9ccc45a63996a967c72983de30"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a14171a9ccc45a63996a967c72983de30">gridfire::AdaptiveEngineView::m_config</a></div><div class="ttdeci">Config &amp; m_config</div><div class="ttdoc">A reference to the singleton Config instance, used for retrieving configuration parameters.</div><div class="ttdef"><b>Definition</b> <a href="#l00263">engine_adaptive.h:263</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a19fc7e02e216b797aa643fa35e429800"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a19fc7e02e216b797aa643fa35e429800">gridfire::AdaptiveEngineView::m_activeReactions</a></div><div class="ttdeci">reaction::LogicalReactionSet m_activeReactions</div><div class="ttdoc">The set of reactions that are currently active in the network.</div><div class="ttdef"><b>Definition</b> <a href="#l00273">engine_adaptive.h:273</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a21c6e33bbf8c18fd5b5eaabb469054de"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a21c6e33bbf8c18fd5b5eaabb469054de">gridfire::AdaptiveEngineView::m_reactionIndexMap</a></div><div class="ttdeci">std::vector&lt; size_t &gt; m_reactionIndexMap</div><div class="ttdoc">A map from the indices of the active reactions to the indices of the corresponding reactions in the f...</div><div class="ttdef"><b>Definition</b> <a href="#l00278">engine_adaptive.h:278</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a231193a61ba5a31e8eb92b0d4ce69111"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a231193a61ba5a31e8eb92b0d4ce69111">gridfire::AdaptiveEngineView::generateStoichiometryMatrix</a></div><div class="ttdeci">void generateStoichiometryMatrix() override</div><div class="ttdoc">Generates the stoichiometry matrix for the active reactions and species.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00160">engine_adaptive.cpp:160</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a256d14a333f9401039b826cc889761a8"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a256d14a333f9401039b826cc889761a8">gridfire::AdaptiveEngineView::mapCulledToFullSpeciesIndex</a></div><div class="ttdeci">size_t mapCulledToFullSpeciesIndex(size_t culledSpeciesIndex) const</div><div class="ttdoc">Maps a culled species index to a full species index.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00242">engine_adaptive.cpp:242</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a3d9d8e862d1c2f0a8ba460c57f6a7f44"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a3d9d8e862d1c2f0a8ba460c57f6a7f44">gridfire::AdaptiveEngineView::mapFullToCulled</a></div><div class="ttdeci">std::vector&lt; double &gt; mapFullToCulled(const std::vector&lt; double &gt; &amp;full) const</div><div class="ttdoc">Maps a vector of full abundances to a vector of culled abundances.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00233">engine_adaptive.cpp:233</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a42417e96fe9fd623458af109401daf08"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a42417e96fe9fd623458af109401daf08">gridfire::AdaptiveEngineView::cullReactionsByFlow</a></div><div class="ttdeci">std::vector&lt; const reaction::LogicalReaction * &gt; cullReactionsByFlow(const std::vector&lt; ReactionFlow &gt; &amp;allFlows, const std::unordered_set&lt; fourdst::atomic::Species &gt; &amp;reachableSpecies, const std::vector&lt; double &gt; &amp;Y_full, double maxFlow) const</div><div class="ttdoc">Culls reactions from the network based on their flow rates.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00342">engine_adaptive.cpp:342</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a4710d218c8a0fd161e994ecd60b48e58"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a4710d218c8a0fd161e994ecd60b48e58">gridfire::AdaptiveEngineView::getJacobianMatrixEntry</a></div><div class="ttdeci">double getJacobianMatrixEntry(const int i_culled, const int j_culled) const override</div><div class="ttdoc">Gets an entry from the Jacobian matrix for the active species.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00149">engine_adaptive.cpp:149</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a4d38b46be9f25c7afe7ddd2b284253f8"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8">gridfire::AdaptiveEngineView::m_baseEngine</a></div><div class="ttdeci">DynamicEngine &amp; m_baseEngine</div><div class="ttdoc">The underlying engine to which this view delegates calculations.</div><div class="ttdef"><b>Definition</b> <a href="#l00268">engine_adaptive.h:268</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a5eaf7c3a4e28cd3a4f34979b88a80103"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a5eaf7c3a4e28cd3a4f34979b88a80103">gridfire::AdaptiveEngineView::LogManager</a></div><div class="ttdeci">fourdst::logging::LogManager LogManager</div><div class="ttdef"><b>Definition</b> <a href="#l00261">engine_adaptive.h:261</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a5f66204a0ff5b27eed243afddecb0093"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a5f66204a0ff5b27eed243afddecb0093">gridfire::AdaptiveEngineView::m_speciesIndexMap</a></div><div class="ttdeci">std::vector&lt; size_t &gt; m_speciesIndexMap</div><div class="ttdoc">A map from the indices of the active species to the indices of the corresponding species in the full ...</div><div class="ttdef"><b>Definition</b> <a href="#l00276">engine_adaptive.h:276</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a63580db57e0f48f508906a11ccfd465e"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a63580db57e0f48f508906a11ccfd465e">gridfire::AdaptiveEngineView::m_isStale</a></div><div class="ttdeci">bool m_isStale</div><div class="ttdoc">A flag indicating whether the view is stale and needs to be updated.</div><div class="ttdef"><b>Definition</b> <a href="#l00281">engine_adaptive.h:281</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a67b4ea8cad115394bb4a42cc39d696f9"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a67b4ea8cad115394bb4a42cc39d696f9">gridfire::AdaptiveEngineView::getStoichiometryMatrixEntry</a></div><div class="ttdeci">int getStoichiometryMatrixEntry(const int speciesIndex_culled, const int reactionIndex_culled) const override</div><div class="ttdoc">Gets an entry from the stoichiometry matrix for the active species and reactions.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00165">engine_adaptive.cpp:165</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a68695f056b660e91285b7e5a931612e1"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a68695f056b660e91285b7e5a931612e1">gridfire::AdaptiveEngineView::mapCulledToFull</a></div><div class="ttdeci">std::vector&lt; double &gt; mapCulledToFull(const std::vector&lt; double &gt; &amp;culled) const</div><div class="ttdoc">Maps a vector of culled abundances to a vector of full abundances.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00224">engine_adaptive.cpp:224</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a7b276b7210be588263395bdb0497fc6d"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a7b276b7210be588263395bdb0497fc6d">gridfire::AdaptiveEngineView::calculateRHSAndEnergy</a></div><div class="ttdeci">StepDerivatives&lt; double &gt; calculateRHSAndEnergy(const std::vector&lt; double &gt; &amp;Y_culled, const double T9, const double rho) const override</div><div class="ttdoc">Calculates the right-hand side (dY/dt) and energy generation for the active species.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00120">engine_adaptive.cpp:120</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a7d61e73f5158f1574cda3edc90c51f7e"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a7d61e73f5158f1574cda3edc90c51f7e">gridfire::AdaptiveEngineView::update</a></div><div class="ttdeci">void update(const NetIn &amp;netIn) override</div><div class="ttdoc">Updates the active species and reactions based on the current conditions.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00085">engine_adaptive.cpp:85</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a89614f4a48f60c4170a0197f45303e7c"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a89614f4a48f60c4170a0197f45303e7c">gridfire::AdaptiveEngineView::constructReactionIndexMap</a></div><div class="ttdeci">std::vector&lt; size_t &gt; constructReactionIndexMap() const</div><div class="ttdoc">Constructs the reaction index map.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00053">engine_adaptive.cpp:53</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a896d29325b4233e83d9298850b617a2d"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a896d29325b4233e83d9298850b617a2d">gridfire::AdaptiveEngineView::constructSpeciesIndexMap</a></div><div class="ttdeci">std::vector&lt; size_t &gt; constructSpeciesIndexMap() const</div><div class="ttdoc">Constructs the species index map.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00024">engine_adaptive.cpp:24</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a91e742642d8a8d9ec0620779927e5101"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a91e742642d8a8d9ec0620779927e5101">gridfire::AdaptiveEngineView::mapCulledToFullReactionIndex</a></div><div class="ttdeci">size_t mapCulledToFullReactionIndex(size_t culledReactionIndex) const</div><div class="ttdoc">Maps a culled reaction index to a full reaction index.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00251">engine_adaptive.cpp:251</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a93b38d0fdc4647f6f7340172dae17872"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a93b38d0fdc4647f6f7340172dae17872">gridfire::AdaptiveEngineView::getSpeciesTimescales</a></div><div class="ttdeci">std::unordered_map&lt; fourdst::atomic::Species, double &gt; getSpeciesTimescales(const std::vector&lt; double &gt; &amp;Y_culled, double T9, double rho) const override</div><div class="ttdoc">Computes timescales for all active species in the network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00196">engine_adaptive.cpp:196</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_aa79fb382c98461b02a2c30668491e6c5"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#aa79fb382c98461b02a2c30668491e6c5">gridfire::AdaptiveEngineView::finalizeActiveSet</a></div><div class="ttdeci">void finalizeActiveSet(const std::vector&lt; const reaction::LogicalReaction * &gt; &amp;finalReactions)</div><div class="ttdoc">Finalizes the set of active species and reactions.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00385">engine_adaptive.cpp:385</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_aae4ddbef1c4e2202fd236221a4bf376b"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#aae4ddbef1c4e2202fd236221a4bf376b">gridfire::AdaptiveEngineView::setScreeningModel</a></div><div class="ttdeci">void setScreeningModel(screening::ScreeningType model) override</div><div class="ttdoc">Sets the screening model for the base engine.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00216">engine_adaptive.cpp:216</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_abdbaf4b87629efe43ac1255dad424c0c"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#abdbaf4b87629efe43ac1255dad424c0c">gridfire::AdaptiveEngineView::calculateAllReactionFlows</a></div><div class="ttdeci">std::vector&lt; ReactionFlow &gt; calculateAllReactionFlows(const NetIn &amp;netIn, std::vector&lt; double &gt; &amp;out_Y_Full) const</div><div class="ttdoc">Calculates the molar reaction flow rate for all reactions in the full network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00268">engine_adaptive.cpp:268</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_ac5bdbe46f87d38d9f23ece5743dcd193"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#ac5bdbe46f87d38d9f23ece5743dcd193">gridfire::AdaptiveEngineView::m_logger</a></div><div class="ttdeci">quill::Logger * m_logger</div><div class="ttdoc">A pointer to the logger instance, used for logging messages.</div><div class="ttdef"><b>Definition</b> <a href="#l00265">engine_adaptive.h:265</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_ac83a8efe25c0e5b9bf7756ac3a500bb1"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#ac83a8efe25c0e5b9bf7756ac3a500bb1">gridfire::AdaptiveEngineView::getNetworkSpecies</a></div><div class="ttdeci">const std::vector&lt; fourdst::atomic::Species &gt; &amp; getNetworkSpecies() const override</div><div class="ttdoc">Gets the list of active species in the network.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00116">engine_adaptive.cpp:116</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_ac9aab6f60e80a9228b2b19b1b10449ef"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#ac9aab6f60e80a9228b2b19b1b10449ef">gridfire::AdaptiveEngineView::generateJacobianMatrix</a></div><div class="ttdeci">void generateJacobianMatrix(const std::vector&lt; double &gt; &amp;Y_culled, const double T9, const double rho) override</div><div class="ttdoc">Generates the Jacobian matrix for the active species.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00138">engine_adaptive.cpp:138</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_ad599363cdd457e72e2e2937b0222c455"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#ad599363cdd457e72e2e2937b0222c455">gridfire::AdaptiveEngineView::AdaptiveEngineView</a></div><div class="ttdeci">AdaptiveEngineView(DynamicEngine &amp;baseEngine)</div><div class="ttdoc">Constructs an AdaptiveEngineView.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00013">engine_adaptive.cpp:13</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_aedc0dedb51c81e03f253cc409a5d5c40"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#aedc0dedb51c81e03f253cc409a5d5c40">gridfire::AdaptiveEngineView::validateState</a></div><div class="ttdeci">void validateState() const</div><div class="ttdoc">Validates that the AdaptiveEngineView is not stale.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00260">engine_adaptive.cpp:260</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_aee095b30a9dce5fcb5ae2fa1d2aa192c"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#aee095b30a9dce5fcb5ae2fa1d2aa192c">gridfire::AdaptiveEngineView::getBaseEngine</a></div><div class="ttdeci">const DynamicEngine &amp; getBaseEngine() const override</div><div class="ttdoc">Gets the base engine.</div><div class="ttdef"><b>Definition</b> <a href="#l00226">engine_adaptive.h:226</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_af4bc9fc6e4afcd6a53c49ca6e2a95940"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#af4bc9fc6e4afcd6a53c49ca6e2a95940">gridfire::AdaptiveEngineView::m_activeSpecies</a></div><div class="ttdeci">std::vector&lt; fourdst::atomic::Species &gt; m_activeSpecies</div><div class="ttdoc">The set of species that are currently active in the network.</div><div class="ttdef"><b>Definition</b> <a href="#l00271">engine_adaptive.h:271</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_afec39b2faa34ea65c5488dd8e11ba3c3"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3">gridfire::AdaptiveEngineView::Config</a></div><div class="ttdeci">fourdst::config::Config Config</div><div class="ttdef"><b>Definition</b> <a href="#l00260">engine_adaptive.h:260</a></div></div>
<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="engine__abstract_8h_source.html#l00121">engine_abstract.h:121</a></div></div>
<div class="ttc" id="aclassgridfire_1_1_engine_view_html"><div class="ttname"><a href="classgridfire_1_1_engine_view.html">gridfire::EngineView</a></div><div class="ttdoc">Abstract base class for a &quot;view&quot; of a reaction network engine.</div><div class="ttdef"><b>Definition</b> <a href="engine__view__abstract_8h_source.html#l00074">engine_view_abstract.h:74</a></div></div>
<div class="ttc" id="aclassgridfire_1_1reaction_1_1_logical_reaction_html"><div class="ttname"><a href="classgridfire_1_1reaction_1_1_logical_reaction.html">gridfire::reaction::LogicalReaction</a></div><div class="ttdoc">Represents a &quot;logical&quot; reaction that aggregates rates from multiple sources.</div><div class="ttdef"><b>Definition</b> <a href="reaction_8h_source.html#l00308">reaction.h:308</a></div></div>
<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="aengine__abstract_8h_html"><div class="ttname"><a href="engine__abstract_8h.html">engine_abstract.h</a></div><div class="ttdoc">Abstract interfaces for reaction network engines in GridFire.</div></div>
<div class="ttc" id="aengine__view__abstract_8h_html"><div class="ttname"><a href="engine__view__abstract_8h.html">engine_view_abstract.h</a></div><div class="ttdoc">Abstract interfaces for engine &quot;views&quot; in GridFire.</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&lt; LogicalReaction &gt; 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="ttdoc">Enumerates the available plasma screening models.</div><div class="ttdef"><b>Definition</b> <a href="screening__types_8h_source.html#l00015">screening_types.h:15</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="engine__abstract_8h_source.html#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>
<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_adaptive_engine_view_1_1_reaction_flow_html"><div class="ttname"><a href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html">gridfire::AdaptiveEngineView::ReactionFlow</a></div><div class="ttdoc">A struct to hold a reaction and its flow rate.</div><div class="ttdef"><b>Definition</b> <a href="#l00287">engine_adaptive.h:287</a></div></div>
<div class="ttc" id="astructgridfire_1_1_adaptive_engine_view_1_1_reaction_flow_html_a3bb21f20df8115d37108cf3c3be3bc6f"><div class="ttname"><a href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a3bb21f20df8115d37108cf3c3be3bc6f">gridfire::AdaptiveEngineView::ReactionFlow::reactionPtr</a></div><div class="ttdeci">const reaction::LogicalReaction * reactionPtr</div><div class="ttdef"><b>Definition</b> <a href="#l00288">engine_adaptive.h:288</a></div></div>
<div class="ttc" id="astructgridfire_1_1_adaptive_engine_view_1_1_reaction_flow_html_a6c2e9087f6c8af5d89a5f0de7bd4a5b4"><div class="ttname"><a href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a6c2e9087f6c8af5d89a5f0de7bd4a5b4">gridfire::AdaptiveEngineView::ReactionFlow::flowRate</a></div><div class="ttdeci">double flowRate</div><div class="ttdef"><b>Definition</b> <a href="#l00289">engine_adaptive.h:289</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="engine__abstract_8h_source.html#l00053">engine_abstract.h:53</a></div></div>
</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="dir_fb341b7e674a7e4701415d4572cba12f.html">views</a></li><li class="navelem"><a class="el" href="engine__adaptive_8h.html">engine_adaptive.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>