Files
GridFire/docs/html/index.html

266 lines
20 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: GridFire</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('index.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><div class="header">
<div class="headertitle"><div class="title">GridFire </div></div>
</div><!--header-->
<div class="contents">
<div class="textblock"><p><a class="anchor" id="md_docs_2static_2mainpage"></a></p>
<p><img src="../../assets/logo/GridFire.png" alt="GridFire Logo" class="inline"/></p>
<h1><a class="anchor" id="autotoc_md0"></a>
GridFire: A Nuclear Reaction Network Library</h1>
<h2><a class="anchor" id="autotoc_md1"></a>
Overview</h2>
<p>GridFire is a C++ library for simulating nuclear reaction networks. It is designed to be flexible, performant, and easy to integrate into larger astrophysical simulations. GridFire provides a modular framework for defining, manipulating, and solving complex reaction networks, with a focus on providing different "engines" to suit various computational needs.</p>
<p>This documentation provides a comprehensive overview of the GridFire library, its architecture, and its components.</p>
<h2><a class="anchor" id="autotoc_md2"></a>
Features</h2>
<ul>
<li><b>Modular Design:</b> The library is organized into distinct components for reactions, engines, and solvers, allowing for easy extension and customization.</li>
<li><b>Multiple Engines:</b> GridFire offers several reaction network engines, including:<ul>
<li>A <code>GraphEngine</code> for solving the full reaction network.</li>
<li>An <code>AdaptiveEngine</code> that dynamically adjusts the network size at runtime for improved performance.</li>
<li>A simplified <code>Approx8Engine</code> for testing and for scenarios where a full network is not required.</li>
</ul>
</li>
<li><b>REACLIB Integration:</b> GridFire can parse and utilize data from the REACLIB database, allowing for the use of standard, up-to-date reaction rates.</li>
<li><b>High Performance:</b> The library is designed for performance, using modern C++ features and high-performance libraries like Eigen, Boost, and CppAD.</li>
</ul>
<h2><a class="anchor" id="autotoc_md3"></a>
Directory Structure</h2>
<p>The GridFire project is organized into the following main directories:</p>
<ul>
<li><code>src/</code>: Contains the source code for the GridFire library.<ul>
<li><code>include/gridfire/</code>: Public header files for the library.<ul>
<li><code>reaction/</code>: Header files related to defining and managing nuclear reactions.</li>
<li><code>engine/</code>: Header files for the various reaction network engines.</li>
<li><code>solver/</code>: Header files for the numerical solvers.</li>
</ul>
</li>
<li><code>lib/</code>: Source code for the implementation of the library components.</li>
</ul>
</li>
<li><code>subprojects/</code>: Contains external libraries and dependencies.</li>
<li><code>tests/</code>: Unit tests for the GridFire library.</li>
<li><code>docs/</code>: Directory for documentation, including this Doxygen homepage.</li>
</ul>
<h2><a class="anchor" id="autotoc_md4"></a>
Core Components</h2>
<p>The GridFire library is built around three main components: <b>Reactions</b>, <b>Engines</b>, and <b>Solvers</b>.</p>
<h3><a class="anchor" id="autotoc_md5"></a>
Reactions</h3>
<p>The <code>reaction</code> component is responsible for defining and managing nuclear reactions. The key classes are:</p>
<ul>
<li><code><a class="el" href="classgridfire_1_1reaction_1_1_reaction.html" title="Represents a single nuclear reaction from a specific data source.">gridfire::reaction::Reaction</a></code>: Represents a single nuclear reaction, including its reactants, products, and rate coefficients.</li>
<li><code><a class="el" href="namespacegridfire_1_1reaction.html#ad838ce3fb6cc02c3fd90b924a0dd91b1" title="A set of reactions, typically from a single source like REACLIB.">gridfire::reaction::ReactionSet</a></code>: A collection of <code>Reaction</code> objects, representing a full reaction network.</li>
<li><code><a class="el" href="classgridfire_1_1reaction_1_1_logical_reaction.html" title="Represents a &quot;logical&quot; reaction that aggregates rates from multiple sources.">gridfire::reaction::LogicalReaction</a></code>: An abstraction that can represent a single reaction or a combination of multiple reactions (e.g., a forward and reverse reaction).</li>
<li><code><a class="el" href="namespacegridfire_1_1reaction.html#aa86f08712565f278adacc7cd2361eb31" title="A set of logical reactions.">gridfire::reaction::LogicalReactionSet</a></code>: A collection of <code>LogicalReaction</code> objects.</li>
</ul>
<p>GridFire can load reaction data from the REACLIB database via the <code>build_reaclib_nuclear_network</code> function, which is implemented in <code><a class="el" href="reaclib_8cpp.html">src/network/lib/reaction/reaclib.cpp</a></code>.</p>
<h3><a class="anchor" id="autotoc_md6"></a>
Engines</h3>
<p>Engines are responsible for the low-level details of solving the reaction network equations. GridFire provides several engines, each with its own trade-offs in terms of accuracy and performance:</p>
<ul>
<li><b><code><a class="el" href="classgridfire_1_1_graph_engine.html" title="A reaction network engine that uses a graph-based representation.">gridfire::GraphEngine</a></code></b>: A full-network solver that represents the reaction network as a graph. It uses the CppAD library for automatic differentiation to compute the Jacobian matrix, which is essential for implicit solvers.</li>
<li><b><code>gridfire::AdaptiveEngine</code></b>: A dynamic engine that adapts the size of the reaction network at runtime. It can add or remove reactions and species based on their importance, which can significantly improve performance in simulations with changing conditions.</li>
<li><b><code>gridfire::Approx8Engine</code></b>: A simplified, 8-isotope network for hydrogen and helium burning. This engine is useful for testing and for simulations where a full network is not necessary. It includes the pp-chain and the CNO cycle.</li>
</ul>
<h3><a class="anchor" id="autotoc_md7"></a>
Solvers</h3>
<p>The <code>solver</code> component provides the numerical algorithms for solving the system of ordinary differential equations (ODEs) that describe the evolution of the species abundances. GridFire uses the <code>boost::numeric::odeint</code> library for its ODE solvers.</p>
<p>The main solver class is <code><a class="el" href="classgridfire_1_1solver_1_1_q_s_e_network_solver.html" title="A network solver that uses a Quasi-Steady-State Equilibrium (QSE) approach.">gridfire::solver::QSENetworkSolver</a></code>, which implements a Quasi-Steady-State (QSE) solver.</p>
<h2><a class="anchor" id="autotoc_md8"></a>
Dependencies</h2>
<p>GridFire relies on the following external libraries:</p>
<ul>
<li><b><a href="https://eigen.tuxfamily.org/">Eigen</a></b>: A C++ template library for linear algebra.</li>
<li><b><a href="https://www.boost.org/">Boost</a></b>: A collection of peer-reviewed, high-quality C++ libraries. GridFire uses Boost for sparse matrices and numerical integration.</li>
<li><b><a href="https://cppad.readthedocs.io/en/latest/">CppAD</a></b>: A C++ template library for automatic differentiation.</li>
<li><b><a href="https://github.com/odygrd/quill">Quill</a></b>: A high-performance, asynchronous logging library.</li>
<li><b><a href="https://github.com/jbeder/yaml-cpp">yaml-cpp</a></b>: A YAML parser and emitter for C++.</li>
<li><b><a href="https://github.com/google/googletest">GoogleTest</a></b>: A unit testing framework for C++.</li>
<li><b><a href="https://github.com/Cyan4973/xxHash">xxHash</a></b>: An extremely fast non-cryptographic hash algorithm.</li>
</ul>
<p>Aside from boost, all dependencies are automatically downloaded and built by the Meson build system. Boost must be installed on your system beforehand.</p>
<h2><a class="anchor" id="autotoc_md9"></a>
Building the Project</h2>
<p>GridFire uses the <a href="https://mesonbuild.com/">Meson</a> build system. To build the project, you will need to have Meson and a C++ compiler (like g++ or clang++) installed.</p>
<p>You can configure the build using the options listed in <code>meson_options.txt</code>.</p>
<div class="fragment"><div class="line">meson setup builddir --prefix=/path/to/install --buildtype=release</div>
<div class="line">meson compile -C builddir</div>
<div class="line">meson install -C builddir</div>
</div><!-- fragment --><h2><a class="anchor" id="autotoc_md10"></a>
How to Use</h2>
<p>To use the GridFire library in your own project, you will need to link against the compiled library and include the necessary header files.</p>
<p>Here is a simple example of how to use the library:</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="network_8h.html">gridfire/network.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="solver_8h.html">gridfire/solver/solver.h</a>&quot;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> main() {</div>
<div class="line"> <span class="comment">// 1. Create a Composition object with the initial abundances</span></div>
<div class="line"> fourdst::composition::Composition comp;</div>
<div class="line"> <span class="comment">// ... set initial abundances ...</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 2. Create a reaction network engine</span></div>
<div class="line"> <a class="code hl_class" href="classgridfire_1_1_graph_engine.html">gridfire::GraphEngine</a> engine(comp);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 3. Create a solver</span></div>
<div class="line"> <a class="code hl_class" href="classgridfire_1_1solver_1_1_q_s_e_network_solver.html">gridfire::solver::QSENetworkSolver</a> solver(engine);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 4. Set the thermodynamic conditions</span></div>
<div class="line"> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">gridfire::NetIn</a> netIn;</div>
<div class="line"> netIn.<a class="code hl_variable" href="structgridfire_1_1_net_in.html#a5be0f5195a5cd1dd177b9fc5ab83a7be">temperature</a> = 1.0e8; <span class="comment">// K</span></div>
<div class="line"> netIn.<a class="code hl_variable" href="structgridfire_1_1_net_in.html#a06f0dff9f8927b7cf2da3004c8fa1577">density</a> = 1.0e4; <span class="comment">// g/cm^3</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 5. Solve the network</span></div>
<div class="line"> <a class="code hl_struct" href="structgridfire_1_1_net_out.html">gridfire::NetOut</a> netOut = solver.evaluate(netIn);</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="aclassgridfire_1_1_graph_engine_html"><div class="ttname"><a href="classgridfire_1_1_graph_engine.html">gridfire::GraphEngine</a></div><div class="ttdoc">A reaction network engine that uses a graph-based representation.</div><div class="ttdef"><b>Definition</b> <a href="engine__graph_8h_source.html#l00089">engine_graph.h:89</a></div></div>
<div class="ttc" id="aclassgridfire_1_1solver_1_1_q_s_e_network_solver_html"><div class="ttname"><a href="classgridfire_1_1solver_1_1_q_s_e_network_solver.html">gridfire::solver::QSENetworkSolver</a></div><div class="ttdoc">A network solver that uses a Quasi-Steady-State Equilibrium (QSE) approach.</div><div class="ttdef"><b>Definition</b> <a href="solver_8h_source.html#l00098">solver.h:98</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="asolver_8h_html"><div class="ttname"><a href="solver_8h.html">solver.h</a></div></div>
<div class="ttc" id="astructgridfire_1_1_net_in_html"><div class="ttname"><a href="structgridfire_1_1_net_in.html">gridfire::NetIn</a></div><div class="ttdef"><b>Definition</b> <a href="network_8h_source.html#l00053">network.h:53</a></div></div>
<div class="ttc" id="astructgridfire_1_1_net_in_html_a06f0dff9f8927b7cf2da3004c8fa1577"><div class="ttname"><a href="structgridfire_1_1_net_in.html#a06f0dff9f8927b7cf2da3004c8fa1577">gridfire::NetIn::density</a></div><div class="ttdeci">double density</div><div class="ttdoc">Density in g/cm^3.</div><div class="ttdef"><b>Definition</b> <a href="network_8h_source.html#l00058">network.h:58</a></div></div>
<div class="ttc" id="astructgridfire_1_1_net_in_html_a5be0f5195a5cd1dd177b9fc5ab83a7be"><div class="ttname"><a href="structgridfire_1_1_net_in.html#a5be0f5195a5cd1dd177b9fc5ab83a7be">gridfire::NetIn::temperature</a></div><div class="ttdeci">double temperature</div><div class="ttdoc">Temperature in Kelvin.</div><div class="ttdef"><b>Definition</b> <a href="network_8h_source.html#l00057">network.h:57</a></div></div>
<div class="ttc" id="astructgridfire_1_1_net_out_html"><div class="ttname"><a href="structgridfire_1_1_net_out.html">gridfire::NetOut</a></div><div class="ttdef"><b>Definition</b> <a href="network_8h_source.html#l00065">network.h:65</a></div></div>
</div><!-- fragment --><p> Linking can be done using the pkg-config tool:</p>
<div class="fragment"><div class="line">pkg-config --cflags --libs gridfire</div>
</div><!-- fragment --><p>Note that you will also need to tell the compiler where to find boost headers.</p>
<p>A more detailed example of how to use an AdaptiveEngine can be found below </p><div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="network_8h.html">gridfire/network.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="solver_8h.html">gridfire/solver/solver.h</a>&quot;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> main() {</div>
<div class="line"> <span class="comment">// 1. Create a Composition object with the initial abundances</span></div>
<div class="line"> fourdst::composition::Composition comp;</div>
<div class="line"> <span class="comment">// ... set initial abundances ...</span></div>
<div class="line"> </div>
<div class="line"> <a class="code hl_class" href="classgridfire_1_1_graph_engine.html">gridfire::GraphEngine</a> baseEngine(comp);</div>
<div class="line"> <a class="code hl_class" href="classgridfire_1_1_adaptive_engine_view.html">gridfire::AdaptiveEngineView</a> adaptiveEngine(comp);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 3. Create a solver</span></div>
<div class="line"> <a class="code hl_class" href="classgridfire_1_1solver_1_1_q_s_e_network_solver.html">gridfire::solver::QSENetworkSolver</a> solver(adaptiveEngine);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 4. Set the thermodynamic conditions</span></div>
<div class="line"> <a class="code hl_struct" href="structgridfire_1_1_net_in.html">gridfire::NetIn</a> netIn;</div>
<div class="line"> netIn.<a class="code hl_variable" href="structgridfire_1_1_net_in.html#a5be0f5195a5cd1dd177b9fc5ab83a7be">temperature</a> = 1.0e8; <span class="comment">// K</span></div>
<div class="line"> netIn.<a class="code hl_variable" href="structgridfire_1_1_net_in.html#a06f0dff9f8927b7cf2da3004c8fa1577">density</a> = 1.0e4; <span class="comment">// g/cm^3</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// 5. Solve the network</span></div>
<div class="line"> <a class="code hl_struct" href="structgridfire_1_1_net_out.html">gridfire::NetOut</a> netOut = solver.evaluate(netIn);</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="aclassgridfire_1_1_adaptive_engine_view_html"><div class="ttname"><a href="classgridfire_1_1_adaptive_engine_view.html">gridfire::AdaptiveEngineView</a></div><div class="ttdoc">An engine view that dynamically adapts the reaction network based on runtime conditions.</div><div class="ttdef"><b>Definition</b> <a href="engine__adaptive_8h_source.html#l00047">engine_adaptive.h:47</a></div></div>
</div><!-- fragment --><p>Note how the adaptive engine is a view of the base engine. This allows the adaptive engine to dynamically adjust the network size at runtime based on the reactions that are active. </p>
</div></div><!-- PageDoc -->
<a href="doxygen_crawl.html"></a>
</div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
<ul>
<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>