332 lines
68 KiB
HTML
332 lines
68 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/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"> 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__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 "<a class="code" href="engine__abstract_8h.html">gridfire/engine/engine_abstract.h</a>"</span></div>
|
|
<div class="line"><a id="l00003" name="l00003"></a><span class="lineno"> 3</span><span class="preprocessor">#include "<a class="code" href="engine__view__abstract_8h.html">gridfire/engine/views/engine_view_abstract.h</a>"</span></div>
|
|
<div class="line"><a id="l00004" name="l00004"></a><span class="lineno"> 4</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="l00005" name="l00005"></a><span class="lineno"> 5</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="l00006" name="l00006"></a><span class="lineno"> 6</span><span class="preprocessor">#include "<a class="code" href="network_8h.html">gridfire/network.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 "fourdst/composition/atomicSpecies.h"</span></div>
|
|
<div class="line"><a id="l00009" name="l00009"></a><span class="lineno"> 9</span><span class="preprocessor">#include "fourdst/config/config.h"</span></div>
|
|
<div class="line"><a id="l00010" name="l00010"></a><span class="lineno"> 10</span><span class="preprocessor">#include "fourdst/logging/logging.h"</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 "<a class="code" href="priming_8h.html">gridfire/engine/procedures/priming.h</a>"</span></div>
|
|
<div class="line"><a id="l00013" name="l00013"></a><span class="lineno"> 13</span><span class="preprocessor">#include "<a class="code" href="construction_8h.html">gridfire/engine/procedures/construction.h</a>"</span></div>
|
|
<div class="line"><a id="l00014" name="l00014"></a><span class="lineno"> 14</span> </div>
|
|
<div class="line"><a id="l00015" name="l00015"></a><span class="lineno"> 15</span><span class="preprocessor">#include "quill/Logger.h"</span></div>
|
|
<div class="line"><a id="l00016" name="l00016"></a><span class="lineno"> 16</span> </div>
|
|
<div class="line"><a id="l00017" name="l00017"></a><span class="lineno"> 17</span><span class="keyword">namespace </span><a class="code hl_namespace" href="namespacegridfire.html">gridfire</a> {</div>
|
|
<div class="foldopen" id="foldopen00050" data-start="{" data-end="};">
|
|
<div class="line"><a id="l00050" name="l00050"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html"> 50</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><DynamicEngine> {</div>
|
|
<div class="line"><a id="l00051" name="l00051"></a><span class="lineno"> 51</span> <span class="keyword">public</span>:</div>
|
|
<div class="line"><a id="l00060" name="l00060"></a><span class="lineno"> 60</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>& baseEngine);</div>
|
|
<div class="line"><a id="l00061" name="l00061"></a><span class="lineno"> 61</span></div>
|
|
<div class="line"><a id="l00079" name="l00079"></a><span class="lineno"> 79</span> fourdst::composition::Composition <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a2a7ecf985a326b4bea43e00cf9ee43dd">update</a>(<span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a> &netIn) <span class="keyword">override</span>;</div>
|
|
<div class="line"><a id="l00080" name="l00080"></a><span class="lineno"> 80</span> </div>
|
|
<div class="line"><a id="l00081" name="l00081"></a><span class="lineno"> 81</span> <span class="keywordtype">bool</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#ad268c9942655e5c9605148fe07718e88">isStale</a>(<span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>& netIn) <span class="keyword">override</span>;</div>
|
|
<div class="line"><a id="l00082" name="l00082"></a><span class="lineno"> 82</span></div>
|
|
<div class="line"><a id="l00087" name="l00087"></a><span class="lineno"> 87</span> [[nodiscard]] <span class="keyword">const</span> std::vector<fourdst::atomic::Species>& <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="l00088" name="l00088"></a><span class="lineno"> 88</span></div>
|
|
<div class="line"><a id="l00105" name="l00105"></a><span class="lineno"> 105</span> [[nodiscard]] std::expected<StepDerivatives<double>, <a class="code hl_struct" href="structgridfire_1_1expectations_1_1_stale_engine_error.html">expectations::StaleEngineError</a>> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#af703ad17ea65ffff4b75bf8ccc00e5d5">calculateRHSAndEnergy</a>(</div>
|
|
<div class="line"><a id="l00106" name="l00106"></a><span class="lineno"> 106</span> <span class="keyword">const</span> std::vector<double> &Y_culled,</div>
|
|
<div class="line"><a id="l00107" name="l00107"></a><span class="lineno"> 107</span> <span class="keyword">const</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00108" name="l00108"></a><span class="lineno"> 108</span> <span class="keyword">const</span> <span class="keywordtype">double</span> rho</div>
|
|
<div class="line"><a id="l00109" name="l00109"></a><span class="lineno"> 109</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00110" name="l00110"></a><span class="lineno"> 110</span></div>
|
|
<div class="line"><a id="l00124" name="l00124"></a><span class="lineno"> 124</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a03fc187d3d306b9058103b9522cbbaeb">generateJacobianMatrix</a>(</div>
|
|
<div class="line"><a id="l00125" name="l00125"></a><span class="lineno"> 125</span> <span class="keyword">const</span> std::vector<double> &Y_dynamic,</div>
|
|
<div class="line"><a id="l00126" name="l00126"></a><span class="lineno"> 126</span> <span class="keyword">const</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00127" name="l00127"></a><span class="lineno"> 127</span> <span class="keyword">const</span> <span class="keywordtype">double</span> rho</div>
|
|
<div class="line"><a id="l00128" name="l00128"></a><span class="lineno"> 128</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00129" name="l00129"></a><span class="lineno"> 129</span></div>
|
|
<div class="line"><a id="l00144" name="l00144"></a><span class="lineno"> 144</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="l00145" name="l00145"></a><span class="lineno"> 145</span> <span class="keyword">const</span> <span class="keywordtype">int</span> i_culled,</div>
|
|
<div class="line"><a id="l00146" name="l00146"></a><span class="lineno"> 146</span> <span class="keyword">const</span> <span class="keywordtype">int</span> j_culled</div>
|
|
<div class="line"><a id="l00147" name="l00147"></a><span class="lineno"> 147</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00148" name="l00148"></a><span class="lineno"> 148</span></div>
|
|
<div class="line"><a id="l00158" name="l00158"></a><span class="lineno"> 158</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="l00159" name="l00159"></a><span class="lineno"> 159</span></div>
|
|
<div class="line"><a id="l00174" name="l00174"></a><span class="lineno"> 174</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="l00175" name="l00175"></a><span class="lineno"> 175</span> <span class="keyword">const</span> <span class="keywordtype">int</span> speciesIndex_culled,</div>
|
|
<div class="line"><a id="l00176" name="l00176"></a><span class="lineno"> 176</span> <span class="keyword">const</span> <span class="keywordtype">int</span> reactionIndex_culled</div>
|
|
<div class="line"><a id="l00177" name="l00177"></a><span class="lineno"> 177</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00178" name="l00178"></a><span class="lineno"> 178</span></div>
|
|
<div class="line"><a id="l00194" name="l00194"></a><span class="lineno"> 194</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="l00195" name="l00195"></a><span class="lineno"> 195</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="l00196" name="l00196"></a><span class="lineno"> 196</span> <span class="keyword">const</span> std::vector<double> &Y_culled,</div>
|
|
<div class="line"><a id="l00197" name="l00197"></a><span class="lineno"> 197</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00198" name="l00198"></a><span class="lineno"> 198</span> <span class="keywordtype">double</span> rho</div>
|
|
<div class="line"><a id="l00199" name="l00199"></a><span class="lineno"> 199</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00200" name="l00200"></a><span class="lineno"> 200</span></div>
|
|
<div class="line"><a id="l00206" name="l00206"></a><span class="lineno"> 206</span> [[nodiscard]] <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_adaptive_engine_view.html#a12cc2f352678fba9688363ba1876ab9c">getNetworkReactions</a>() <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00207" name="l00207"></a><span class="lineno"> 207</span> </div>
|
|
<div class="line"><a id="l00208" name="l00208"></a><span class="lineno"> 208</span> <span class="keywordtype">void</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a7b3a6b3ab0a52f0f84d2b142e74ea672">setNetworkReactions</a>(<span class="keyword">const</span> <a class="code hl_typedef" href="namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31">reaction::LogicalReactionSet</a>& reactions) <span class="keyword">override</span>;</div>
|
|
<div class="line"><a id="l00209" name="l00209"></a><span class="lineno"> 209</span></div>
|
|
<div class="line"><a id="l00223" name="l00223"></a><span class="lineno"> 223</span> [[nodiscard]] std::expected<std::unordered_map<fourdst::atomic::Species, double>, <a class="code hl_struct" href="structgridfire_1_1expectations_1_1_stale_engine_error.html">expectations::StaleEngineError</a>> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a4e856d6d4d2fc220952bbb7e6b2f85d9">getSpeciesTimescales</a>(</div>
|
|
<div class="line"><a id="l00224" name="l00224"></a><span class="lineno"> 224</span> <span class="keyword">const</span> std::vector<double> &Y_culled,</div>
|
|
<div class="line"><a id="l00225" name="l00225"></a><span class="lineno"> 225</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00226" name="l00226"></a><span class="lineno"> 226</span> <span class="keywordtype">double</span> rho</div>
|
|
<div class="line"><a id="l00227" name="l00227"></a><span class="lineno"> 227</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00228" name="l00228"></a><span class="lineno"> 228</span> </div>
|
|
<div class="line"><a id="l00229" name="l00229"></a><span class="lineno"> 229</span> [[nodiscard]] std::expected<std::unordered_map<fourdst::atomic::Species, double>, <a class="code hl_struct" href="structgridfire_1_1expectations_1_1_stale_engine_error.html">expectations::StaleEngineError</a>> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a522e78bce9ff062939572248d57f8cea">getSpeciesDestructionTimescales</a>(</div>
|
|
<div class="line"><a id="l00230" name="l00230"></a><span class="lineno"> 230</span> <span class="keyword">const</span> std::vector<double> &Y,</div>
|
|
<div class="line"><a id="l00231" name="l00231"></a><span class="lineno"> 231</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00232" name="l00232"></a><span class="lineno"> 232</span> <span class="keywordtype">double</span> rho</div>
|
|
<div class="line"><a id="l00233" name="l00233"></a><span class="lineno"> 233</span> ) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00234" name="l00234"></a><span class="lineno"> 234</span></div>
|
|
<div class="line"><a id="l00239" name="l00239"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#aee095b30a9dce5fcb5ae2fa1d2aa192c"> 239</a></span> [[nodiscard]] <span class="keyword">const</span> <a class="code hl_class" href="classgridfire_1_1_dynamic_engine.html">DynamicEngine</a>& <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="l00240" name="l00240"></a><span class="lineno"> 240</span></div>
|
|
<div class="line"><a id="l00256" name="l00256"></a><span class="lineno"> 256</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="l00257" name="l00257"></a><span class="lineno"> 257</span></div>
|
|
<div class="line"><a id="l00271" name="l00271"></a><span class="lineno"> 271</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="l00272" name="l00272"></a><span class="lineno"> 272</span> </div>
|
|
<div class="line"><a id="l00273" name="l00273"></a><span class="lineno"> 273</span> [[nodiscard]] <span class="keywordtype">int</span> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a9055feb245524a5a9549ace935f059ff">getSpeciesIndex</a>(<span class="keyword">const</span> fourdst::atomic::Species &species) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00274" name="l00274"></a><span class="lineno"> 274</span> </div>
|
|
<div class="line"><a id="l00275" name="l00275"></a><span class="lineno"> 275</span> [[nodiscard]] std::vector<double> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a7d0237956bf3ec7230bc51d88e7f8019">mapNetInToMolarAbundanceVector</a>(<span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a> &netIn) <span class="keyword">const override</span>;</div>
|
|
<div class="line"><a id="l00276" name="l00276"></a><span class="lineno"> 276</span> </div>
|
|
<div class="line"><a id="l00277" name="l00277"></a><span class="lineno"> 277</span> [[nodiscard]] <a class="code hl_struct" href="structgridfire_1_1_priming_report.html">PrimingReport</a> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a70005361262bc180d4417b608661e3c3">primeEngine</a>(<span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a> &netIn) <span class="keyword">override</span>;</div>
|
|
<div class="line"><a id="l00278" name="l00278"></a><span class="lineno"> 278</span> <span class="keyword">private</span>:</div>
|
|
<div class="line"><a id="l00279" name="l00279"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3"> 279</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="l00280" name="l00280"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a5eaf7c3a4e28cd3a4f34979b88a80103"> 280</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="l00282" name="l00282"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a14171a9ccc45a63996a967c72983de30"> 282</a></span> <a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#afec39b2faa34ea65c5488dd8e11ba3c3">Config</a>& <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="l00284" name="l00284"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#ac5bdbe46f87d38d9f23ece5743dcd193"> 284</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">"log"</span>);</div>
|
|
<div class="line"><a id="l00285" name="l00285"></a><span class="lineno"> 285</span></div>
|
|
<div class="line"><a id="l00287" name="l00287"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8"> 287</a></span> <a class="code hl_class" href="classgridfire_1_1_dynamic_engine.html">DynamicEngine</a>& <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a4d38b46be9f25c7afe7ddd2b284253f8">m_baseEngine</a>;</div>
|
|
<div class="line"><a id="l00288" name="l00288"></a><span class="lineno"> 288</span></div>
|
|
<div class="line"><a id="l00290" name="l00290"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#af4bc9fc6e4afcd6a53c49ca6e2a95940"> 290</a></span> std::vector<fourdst::atomic::Species> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#af4bc9fc6e4afcd6a53c49ca6e2a95940">m_activeSpecies</a>;</div>
|
|
<div class="line"><a id="l00292" name="l00292"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a19fc7e02e216b797aa643fa35e429800"> 292</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="l00293" name="l00293"></a><span class="lineno"> 293</span></div>
|
|
<div class="line"><a id="l00295" name="l00295"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a5f66204a0ff5b27eed243afddecb0093"> 295</a></span> std::vector<size_t> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a5f66204a0ff5b27eed243afddecb0093">m_speciesIndexMap</a>;</div>
|
|
<div class="line"><a id="l00297" name="l00297"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a21c6e33bbf8c18fd5b5eaabb469054de"> 297</a></span> std::vector<size_t> <a class="code hl_variable" href="classgridfire_1_1_adaptive_engine_view.html#a21c6e33bbf8c18fd5b5eaabb469054de">m_reactionIndexMap</a>;</div>
|
|
<div class="line"><a id="l00298" name="l00298"></a><span class="lineno"> 298</span></div>
|
|
<div class="line"><a id="l00300" name="l00300"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a63580db57e0f48f508906a11ccfd465e"> 300</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="l00301" name="l00301"></a><span class="lineno"> 301</span> </div>
|
|
<div class="line"><a id="l00302" name="l00302"></a><span class="lineno"> 302</span> <span class="keyword">private</span>:</div>
|
|
<div class="foldopen" id="foldopen00306" data-start="{" data-end="};">
|
|
<div class="line"><a id="l00306" name="l00306"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html"> 306</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="l00307" name="l00307"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a3bb21f20df8115d37108cf3c3be3bc6f"> 307</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="l00308" name="l00308"></a><span class="lineno"><a class="line" href="structgridfire_1_1_adaptive_engine_view_1_1_reaction_flow.html#a6c2e9087f6c8af5d89a5f0de7bd4a5b4"> 308</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="l00309" name="l00309"></a><span class="lineno"> 309</span> };</div>
|
|
</div>
|
|
<div class="line"><a id="l00310" name="l00310"></a><span class="lineno"> 310</span> <span class="keyword">private</span>:</div>
|
|
<div class="line"><a id="l00321" name="l00321"></a><span class="lineno"> 321</span> [[nodiscard]] std::vector<size_t> <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="l00322" name="l00322"></a><span class="lineno"> 322</span></div>
|
|
<div class="line"><a id="l00333" name="l00333"></a><span class="lineno"> 333</span> [[nodiscard]] std::vector<size_t> <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="l00334" name="l00334"></a><span class="lineno"> 334</span></div>
|
|
<div class="line"><a id="l00342" name="l00342"></a><span class="lineno"> 342</span> [[nodiscard]] std::vector<double> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a68695f056b660e91285b7e5a931612e1">mapCulledToFull</a>(<span class="keyword">const</span> std::vector<double>& culled) <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="l00351" name="l00351"></a><span class="lineno"> 351</span> [[nodiscard]] std::vector<double> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a3d9d8e862d1c2f0a8ba460c57f6a7f44">mapFullToCulled</a>(<span class="keyword">const</span> std::vector<double>& full) <span class="keyword">const</span>;</div>
|
|
<div class="line"><a id="l00352" name="l00352"></a><span class="lineno"> 352</span></div>
|
|
<div class="line"><a id="l00361" name="l00361"></a><span class="lineno"> 361</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="l00362" name="l00362"></a><span class="lineno"> 362</span></div>
|
|
<div class="line"><a id="l00371" name="l00371"></a><span class="lineno"> 371</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="l00372" name="l00372"></a><span class="lineno"> 372</span></div>
|
|
<div class="line"><a id="l00378" name="l00378"></a><span class="lineno"> 378</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="l00379" name="l00379"></a><span class="lineno"> 379</span></div>
|
|
<div class="line"><a id="l00401" name="l00401"></a><span class="lineno"> 401</span> std::vector<ReactionFlow> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#abdbaf4b87629efe43ac1255dad424c0c">calculateAllReactionFlows</a>(</div>
|
|
<div class="line"><a id="l00402" name="l00402"></a><span class="lineno"> 402</span> <span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>& netIn,</div>
|
|
<div class="line"><a id="l00403" name="l00403"></a><span class="lineno"> 403</span> std::vector<double>& out_Y_Full</div>
|
|
<div class="line"><a id="l00404" name="l00404"></a><span class="lineno"> 404</span> ) <span class="keyword">const</span>;</div>
|
|
<div class="line"><a id="l00422" name="l00422"></a><span class="lineno"> 422</span> [[nodiscard]] std::unordered_set<fourdst::atomic::Species> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a0ed21f7e7c1034fc87b40d4116c4221b">findReachableSpecies</a>(</div>
|
|
<div class="line"><a id="l00423" name="l00423"></a><span class="lineno"> 423</span> <span class="keyword">const</span> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">NetIn</a>& netIn</div>
|
|
<div class="line"><a id="l00424" name="l00424"></a><span class="lineno"> 424</span> ) <span class="keyword">const</span>;</div>
|
|
<div class="line"><a id="l00445" name="l00445"></a><span class="lineno"> 445</span> [[nodiscard]] std::vector<const reaction::LogicalReaction*> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a42417e96fe9fd623458af109401daf08">cullReactionsByFlow</a>(</div>
|
|
<div class="line"><a id="l00446" name="l00446"></a><span class="lineno"> 446</span> <span class="keyword">const</span> std::vector<ReactionFlow>& allFlows,</div>
|
|
<div class="line"><a id="l00447" name="l00447"></a><span class="lineno"> 447</span> <span class="keyword">const</span> std::unordered_set<fourdst::atomic::Species>& reachableSpecies,</div>
|
|
<div class="line"><a id="l00448" name="l00448"></a><span class="lineno"> 448</span> <span class="keyword">const</span> std::vector<double>& Y_full,</div>
|
|
<div class="line"><a id="l00449" name="l00449"></a><span class="lineno"> 449</span> <span class="keywordtype">double</span> maxFlow</div>
|
|
<div class="line"><a id="l00450" name="l00450"></a><span class="lineno"> 450</span> ) <span class="keyword">const</span>;</div>
|
|
<div class="line"><a id="l00451" name="l00451"></a><span class="lineno"> 451</span> </div>
|
|
<div class="line"><a id="l00452" name="l00452"></a><span class="lineno"><a class="line" href="classgridfire_1_1_adaptive_engine_view.html#a4ff60b5214ec0bdaf683feb6615573a5"> 452</a></span> <span class="keyword">typedef</span> std::pair<std::unordered_set<const reaction::LogicalReaction*>, std::unordered_set<fourdst::atomic::Species>> <a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#a4ff60b5214ec0bdaf683feb6615573a5">RescueSet</a>;</div>
|
|
<div class="line"><a id="l00453" name="l00453"></a><span class="lineno"> 453</span> [[nodiscard]] <a class="code hl_typedef" href="classgridfire_1_1_adaptive_engine_view.html#a4ff60b5214ec0bdaf683feb6615573a5">RescueSet</a> <a class="code hl_function" href="classgridfire_1_1_adaptive_engine_view.html#a823c665ba89452aa2b3d5422fa5d313f">rescueEdgeSpeciesDestructionChannel</a>(</div>
|
|
<div class="line"><a id="l00454" name="l00454"></a><span class="lineno"> 454</span> <span class="keyword">const</span> std::vector<double>& Y_full,</div>
|
|
<div class="line"><a id="l00455" name="l00455"></a><span class="lineno"> 455</span> <span class="keyword">const</span> <span class="keywordtype">double</span> T9,</div>
|
|
<div class="line"><a id="l00456" name="l00456"></a><span class="lineno"> 456</span> <span class="keyword">const</span> <span class="keywordtype">double</span> rho,</div>
|
|
<div class="line"><a id="l00457" name="l00457"></a><span class="lineno"> 457</span> <span class="keyword">const</span> std::vector<fourdst::atomic::Species>& activeSpecies,</div>
|
|
<div class="line"><a id="l00458" name="l00458"></a><span class="lineno"> 458</span> <span class="keyword">const</span> <a class="code hl_typedef" href="namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31">reaction::LogicalReactionSet</a>& activeReactions</div>
|
|
<div class="line"><a id="l00459" name="l00459"></a><span class="lineno"> 459</span> ) <span class="keyword">const</span>;</div>
|
|
<div class="line"><a id="l00475" name="l00475"></a><span class="lineno"> 475</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="l00476" name="l00476"></a><span class="lineno"> 476</span> <span class="keyword">const</span> std::vector<const reaction::LogicalReaction*>& finalReactions</div>
|
|
<div class="line"><a id="l00477" name="l00477"></a><span class="lineno"> 477</span> );</div>
|
|
<div class="line"><a id="l00478" name="l00478"></a><span class="lineno"> 478</span> };</div>
|
|
</div>
|
|
<div class="line"><a id="l00479" name="l00479"></a><span class="lineno"> 479</span>}</div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a03fc187d3d306b9058103b9522cbbaeb"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a03fc187d3d306b9058103b9522cbbaeb">gridfire::AdaptiveEngineView::generateJacobianMatrix</a></div><div class="ttdeci">void generateJacobianMatrix(const std::vector< double > &Y_dynamic, const double T9, const double rho) const 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#l00175">engine_adaptive.cpp:175</a></div></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 &reaction, const std::vector< double > &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#l00212">engine_adaptive.cpp:212</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#l00296">engine_adaptive.cpp:296</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< fourdst::atomic::Species > findReachableSpecies(const NetIn &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#l00399">engine_adaptive.cpp:399</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 & 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#l00229">engine_adaptive.cpp:229</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 & 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="#l00282">engine_adaptive.h:282</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="#l00292">engine_adaptive.h:292</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< size_t > 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="#l00297">engine_adaptive.h:297</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#l00197">engine_adaptive.cpp:197</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#l00341">engine_adaptive.cpp:341</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a2a7ecf985a326b4bea43e00cf9ee43dd"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a2a7ecf985a326b4bea43e00cf9ee43dd">gridfire::AdaptiveEngineView::update</a></div><div class="ttdeci">fourdst::composition::Composition update(const NetIn &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#l00088">engine_adaptive.cpp:88</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< double > mapFullToCulled(const std::vector< double > &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#l00332">engine_adaptive.cpp:332</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< const reaction::LogicalReaction * > cullReactionsByFlow(const std::vector< ReactionFlow > &allFlows, const std::unordered_set< fourdst::atomic::Species > &reachableSpecies, const std::vector< double > &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#l00442">engine_adaptive.cpp:442</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#l00186">engine_adaptive.cpp:186</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 & m_baseEngine</div><div class="ttdoc">The underlying engine to which this view delegates calculations.</div><div class="ttdef"><b>Definition</b> <a href="#l00287">engine_adaptive.h:287</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a4e856d6d4d2fc220952bbb7e6b2f85d9"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a4e856d6d4d2fc220952bbb7e6b2f85d9">gridfire::AdaptiveEngineView::getSpeciesTimescales</a></div><div class="ttdeci">std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesTimescales(const std::vector< double > &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#l00238">engine_adaptive.cpp:238</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a4ff60b5214ec0bdaf683feb6615573a5"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a4ff60b5214ec0bdaf683feb6615573a5">gridfire::AdaptiveEngineView::RescueSet</a></div><div class="ttdeci">std::pair< std::unordered_set< const reaction::LogicalReaction * >, std::unordered_set< fourdst::atomic::Species > > RescueSet</div><div class="ttdef"><b>Definition</b> <a href="#l00452">engine_adaptive.h:452</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a522e78bce9ff062939572248d57f8cea"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a522e78bce9ff062939572248d57f8cea">gridfire::AdaptiveEngineView::getSpeciesDestructionTimescales</a></div><div class="ttdeci">std::expected< std::unordered_map< fourdst::atomic::Species, double >, expectations::StaleEngineError > getSpeciesDestructionTimescales(const std::vector< double > &Y, double T9, double rho) const override</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00266">engine_adaptive.cpp:266</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="#l00280">engine_adaptive.h:280</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< size_t > 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="#l00295">engine_adaptive.h:295</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="#l00300">engine_adaptive.h:300</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#l00202">engine_adaptive.cpp:202</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< double > mapCulledToFull(const std::vector< double > &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#l00323">engine_adaptive.cpp:323</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a70005361262bc180d4417b608661e3c3"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a70005361262bc180d4417b608661e3c3">gridfire::AdaptiveEngineView::primeEngine</a></div><div class="ttdeci">PrimingReport primeEngine(const NetIn &netIn) override</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00308">engine_adaptive.cpp:308</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a7b3a6b3ab0a52f0f84d2b142e74ea672"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a7b3a6b3ab0a52f0f84d2b142e74ea672">gridfire::AdaptiveEngineView::setNetworkReactions</a></div><div class="ttdeci">void setNetworkReactions(const reaction::LogicalReactionSet &reactions) override</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_a7d0237956bf3ec7230bc51d88e7f8019"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a7d0237956bf3ec7230bc51d88e7f8019">gridfire::AdaptiveEngineView::mapNetInToMolarAbundanceVector</a></div><div class="ttdeci">std::vector< double > mapNetInToMolarAbundanceVector(const NetIn &netIn) const override</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00300">engine_adaptive.cpp:300</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a823c665ba89452aa2b3d5422fa5d313f"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a823c665ba89452aa2b3d5422fa5d313f">gridfire::AdaptiveEngineView::rescueEdgeSpeciesDestructionChannel</a></div><div class="ttdeci">RescueSet rescueEdgeSpeciesDestructionChannel(const std::vector< double > &Y_full, const double T9, const double rho, const std::vector< fourdst::atomic::Species > &activeSpecies, const reaction::LogicalReactionSet &activeReactions) const</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00485">engine_adaptive.cpp:485</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< size_t > 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#l00056">engine_adaptive.cpp:56</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< size_t > 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#l00027">engine_adaptive.cpp:27</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_a9055feb245524a5a9549ace935f059ff"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#a9055feb245524a5a9549ace935f059ff">gridfire::AdaptiveEngineView::getSpeciesIndex</a></div><div class="ttdeci">int getSpeciesIndex(const fourdst::atomic::Species &species) const override</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00312">engine_adaptive.cpp:312</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#l00350">engine_adaptive.cpp:350</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< const reaction::LogicalReaction * > &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#l00641">engine_adaptive.cpp:641</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#l00292">engine_adaptive.cpp:292</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< ReactionFlow > calculateAllReactionFlows(const NetIn &netIn, std::vector< double > &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#l00368">engine_adaptive.cpp:368</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="#l00284">engine_adaptive.h:284</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< fourdst::atomic::Species > & 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#l00148">engine_adaptive.cpp:148</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_ad268c9942655e5c9605148fe07718e88"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#ad268c9942655e5c9605148fe07718e88">gridfire::AdaptiveEngineView::isStale</a></div><div class="ttdeci">bool isStale(const NetIn &netIn) override</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00144">engine_adaptive.cpp:144</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 &baseEngine)</div><div class="ttdoc">Constructs an AdaptiveEngineView.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8cpp_source.html#l00016">engine_adaptive.cpp:16</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#l00359">engine_adaptive.cpp:359</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 & getBaseEngine() const override</div><div class="ttdoc">Gets the base engine.</div><div class="ttdef"><b>Definition</b> <a href="#l00239">engine_adaptive.h:239</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< fourdst::atomic::Species > 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="#l00290">engine_adaptive.h:290</a></div></div>
|
|
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html_af703ad17ea65ffff4b75bf8ccc00e5d5"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html#af703ad17ea65ffff4b75bf8ccc00e5d5">gridfire::AdaptiveEngineView::calculateRHSAndEnergy</a></div><div class="ttdeci">std::expected< StepDerivatives< double >, expectations::StaleEngineError > calculateRHSAndEnergy(const std::vector< double > &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#l00152">engine_adaptive.cpp:152</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="#l00279">engine_adaptive.h:279</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#l00130">engine_abstract.h:130</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 "view" 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 "logical" reaction that aggregates rates from multiple sources.</div><div class="ttdef"><b>Definition</b> <a href="reaction_8h_source.html#l00310">reaction.h:310</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="aconstruction_8h_html"><div class="ttname"><a href="construction_8h.html">construction.h</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 "views" 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< LogicalReaction > LogicalReactionSet</div><div class="ttdoc">A set of logical reactions.</div><div class="ttdef"><b>Definition</b> <a href="reaction_8h_source.html#l00563">reaction.h:563</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#l00031">engine_abstract.h:31</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="apriming_8h_html"><div class="ttname"><a href="priming_8h.html">priming.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="#l00306">engine_adaptive.h:306</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="#l00307">engine_adaptive.h:307</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="#l00308">engine_adaptive.h:308</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_priming_report_html"><div class="ttname"><a href="structgridfire_1_1_priming_report.html">gridfire::PrimingReport</a></div><div class="ttdoc">Captures the result of a network priming operation.</div><div class="ttdef"><b>Definition</b> <a href="reporting_8h_source.html#l00067">reporting.h:67</a></div></div>
|
|
<div class="ttc" id="astructgridfire_1_1expectations_1_1_stale_engine_error_html"><div class="ttname"><a href="structgridfire_1_1expectations_1_1_stale_engine_error.html">gridfire::expectations::StaleEngineError</a></div><div class="ttdef"><b>Definition</b> <a href="expected__engine_8h_source.html#l00044">expected_engine.h:44</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_b0856f6b0d80ccb263b2f415c91f9e17.html">include</a></li><li class="navelem"><a class="el" href="dir_3626e0c0e3c5d7812d6b277dfa4ec364.html">gridfire</a></li><li class="navelem"><a class="el" href="dir_aff155d61c3b73b9ab7dcdc908c4d49e.html">engine</a></li><li class="navelem"><a class="el" href="dir_d5492b42d970deba31f48df1b35a6c47.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>
|