Files
SERiF/docs/html/index.html

397 lines
33 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>SERiF: 4DSSE: A 4D Stellar Structure and Evolution Code</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">SERiF<span id="projectnumber">&#160;0.0.1a</span>
</div>
<div id="projectbrief">3+1D Stellar Structure and Evolution</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">4DSSE: A 4D Stellar Structure and Evolution Code </div></div>
</div><!--header-->
<div class="contents">
<div class="textblock"><h1><a class="anchor" id="intro_sec"></a>
Introduction</h1>
<p>Welcome to the documentation for 4DSSE (4D Stellar Structure and Evolution), a new code designed for simulating stellar phenomena in three spatial dimensions plus time. This project is currently under active development.</p>
<p>The primary goal of 4DSSE is to provide a flexible and extensible framework for advanced stellar modeling, incorporating modern numerical techniques and physics modules.</p>
<h1><a class="anchor" id="build_sec"></a>
Building 4DSSE</h1>
<p>The project uses Meson as its build system. MFEM is a core dependency and will be automatically downloaded and built if not found.</p>
<p><b>Prerequisites:</b></p><ul>
<li>A C++ compiler (supporting C++17 or later)</li>
<li>Meson (Install via pip: <code>pip install meson</code>)</li>
<li>Python 3</li>
</ul>
<p><b>Build Steps:</b></p>
<ol type="1">
<li><b>Clone the repository (if you haven't already):</b> <code>bash git clone &lt;repository-url&gt; cd 4dsse </code></li>
<li><b>Using the <code>mk</code> script (recommended for ease of use):</b><ul>
<li>To build the project: <code>bash ./mk </code></li>
<li>To build without tests: <code>bash ./mk --noTest </code></li>
</ul>
</li>
<li><b>Using the <code>4DSSEConsole.sh</code> script:</b> This script provides a simple interface for building and debugging. <code>bash ./4DSSEConsole.sh </code> Follow the on-screen prompts.</li>
<li><b>Manual Meson build:</b><ul>
<li>Setup the build directory (e.g., <code>build</code>): <code>bash meson setup build </code></li>
<li>Compile the project: <code>bash meson compile -C build </code></li>
<li>To run tests (if built with tests enabled): <code>bash meson test -C build </code></li>
</ul>
</li>
</ol>
<p>The compiled executables and libraries will typically be found in the <code>build</code> directory.</p>
<h1><a class="anchor" id="usage_sec"></a>
High-Level Usage Examples</h1>
<p>Below are some high-level examples of how to use key components of 4DSSE.</p>
<h2><a class="anchor" id="usage_polysolver"></a>
Solving for a Polytrope</h2>
<p>The <code>PolySolver</code> class handles the setup and solution of the Lane-Emden equation for polytropic stellar models.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="poly_solver_8h.html">polySolver.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="config_8h.html">config.h</a>&quot;</span> <span class="comment">// For global configuration</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="probe_8h.html">probe.h</a>&quot;</span> <span class="comment">// For logging</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> <span class="comment">// Initialize configuration and logging</span></div>
<div class="line"> Config::getInstance().loadConfig(<span class="stringliteral">&quot;path/to/your/config.yaml&quot;</span>);</div>
<div class="line"> Probe::LogManager::getInstance().getLogger(<span class="stringliteral">&quot;main_log&quot;</span>); <span class="comment">// Initialize a logger</span></div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">try</span> {</div>
<div class="line"> <span class="comment">// Create a PolySolver for polytropic index n=1.5 and FE order 2</span></div>
<div class="line"> PolySolver solver(1.5, 2);</div>
<div class="line"> solver.solve(); <span class="comment">// Solve the system</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Access the solution</span></div>
<div class="line"> mfem::GridFunction&amp; theta_solution = solver.getSolution();</div>
<div class="line"> <span class="comment">// ... process or visualize theta_solution ...</span></div>
<div class="line"> </div>
<div class="line"> } <span class="keywordflow">catch</span> (<span class="keyword">const</span> std::exception&amp; e) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;An error occurred: &quot;</span> &lt;&lt; e.what() &lt;&lt; std::endl;</div>
<div class="line"> <span class="keywordflow">return</span> 1;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="acomp_8cpp_html_ac4c0f8a8146b128f1b8f920e3a9c3b1e"><div class="ttname"><a href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a></div><div class="ttdeci">int main(int argv, char *argc[])</div><div class="ttdef"><b>Definition</b> <a href="comp_8cpp_source.html#l00005">comp.cpp:5</a></div></div>
<div class="ttc" id="aconfig_8h_html"><div class="ttname"><a href="config_8h.html">config.h</a></div></div>
<div class="ttc" id="apoly_solver_8h_html"><div class="ttname"><a href="poly_solver_8h.html">polySolver.h</a></div></div>
<div class="ttc" id="aprobe_8h_html"><div class="ttname"><a href="probe_8h.html">probe.h</a></div></div>
</div><!-- fragment --><h2><a class="anchor" id="usage_composition"></a>
Managing Chemical Compositions</h2>
<p>The <code>Composition</code> class allows for defining and managing chemical compositions.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="composition_8h.html">composition.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;vector&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> <span class="keywordflow">try</span> {</div>
<div class="line"> <span class="comment">// Define symbols and their mass fractions</span></div>
<div class="line"> std::vector&lt;std::string&gt; symbols = {<span class="stringliteral">&quot;H-1&quot;</span>, <span class="stringliteral">&quot;He-4&quot;</span>}; <span class="comment">// Use specific isotopes</span></div>
<div class="line"> std::vector&lt;double&gt; mass_fractions = {0.75, 0.25};</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Create and finalize the composition</span></div>
<div class="line"> composition::Composition comp(symbols, mass_fractions);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get mass fraction of a specific element</span></div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Mass fraction of H-1: &quot;</span> &lt;&lt; comp.getMassFraction(<span class="stringliteral">&quot;H-1&quot;</span>) &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get global properties</span></div>
<div class="line"> <span class="keyword">auto</span> global_props = comp.getComposition().second;</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Mean particle mass: &quot;</span> &lt;&lt; global_props.meanParticleMass &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> } <span class="keywordflow">catch</span> (<span class="keyword">const</span> std::exception&amp; e) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Composition error: &quot;</span> &lt;&lt; e.what() &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="acomposition_8h_html"><div class="ttname"><a href="composition_8h.html">composition.h</a></div></div>
</div><!-- fragment --><h2><a class="anchor" id="usage_network"></a>
Nuclear Reaction Networks</h2>
<p>The <code>Network</code> and <code>Approx8Network</code> classes provide interfaces for nuclear reaction network calculations.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="network_8h.html">network.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="approx8_8h.html">approx8.h</a>&quot;</span> <span class="comment">// Specific network implementation</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;vector&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> nnApprox8::Approx8Network approx8_net; <span class="comment">// Using the Approx8 network</span></div>
<div class="line"> approx8_net.setStiff(<span class="keyword">true</span>); <span class="comment">// Example: use stiff solver</span></div>
<div class="line"> </div>
<div class="line"> nuclearNetwork::NetIn input;</div>
<div class="line"> input.composition = {0.7, 0.0, 0.28, 0.01, 0.005, 0.004, 0.0005, 0.0005}; <span class="comment">// H1, He3, He4, C12, N14, O16, Ne20, Mg24</span></div>
<div class="line"> input.temperature = 1.5e7; <span class="comment">// K</span></div>
<div class="line"> input.density = 150.0; <span class="comment">// g/cm^3</span></div>
<div class="line"> input.tmax = 1.0e10; <span class="comment">// s</span></div>
<div class="line"> input.dt0 = 1.0e6; <span class="comment">// s</span></div>
<div class="line"> <span class="comment">// input.energy can also be set if needed</span></div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">try</span> {</div>
<div class="line"> nuclearNetwork::NetOut output = approx8_net.evaluate(input);</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Number of steps: &quot;</span> &lt;&lt; output.num_steps &lt;&lt; std::endl;</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Final H-1 mass fraction (approx): &quot;</span> &lt;&lt; output.composition[nnApprox8::Net::ih1] &lt;&lt; std::endl;</div>
<div class="line"> } <span class="keywordflow">catch</span> (<span class="keyword">const</span> std::exception&amp; e) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Network evaluation error: &quot;</span> &lt;&lt; e.what() &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="aapprox8_8h_html"><div class="ttname"><a href="approx8_8h.html">approx8.h</a></div><div class="ttdoc">Header file for the Approx8 nuclear reaction network.</div></div>
<div class="ttc" id="anetwork_8h_html"><div class="ttname"><a href="network_8h.html">network.h</a></div></div>
</div><!-- fragment --><h2><a class="anchor" id="usage_constants"></a>
Accessing Physical Constants</h2>
<p>The <code>Constants</code> singleton provides access to a database of physical constants.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="const_8h.html">const.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> Constants&amp; consts = Constants::getInstance();</div>
<div class="line"> <span class="keywordflow">if</span> (consts.isLoaded()) {</div>
<div class="line"> Constant G = consts.get(<span class="stringliteral">&quot;Gravitational constant&quot;</span>);</div>
<div class="line"> std::cout &lt;&lt; G.name &lt;&lt; <span class="stringliteral">&quot;: &quot;</span> &lt;&lt; G.value &lt;&lt; <span class="stringliteral">&quot; &quot;</span> &lt;&lt; G.unit &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> Constant c = consts[<span class="stringliteral">&quot;Speed of light in vacuum&quot;</span>]; <span class="comment">// Can also use operator[]</span></div>
<div class="line"> std::cout &lt;&lt; c.name &lt;&lt; <span class="stringliteral">&quot;: &quot;</span> &lt;&lt; c.value &lt;&lt; <span class="stringliteral">&quot; &quot;</span> &lt;&lt; c.unit &lt;&lt; std::endl;</div>
<div class="line"> } <span class="keywordflow">else</span> {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Failed to load constants.&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="aconst_8h_html"><div class="ttname"><a href="const_8h.html">const.h</a></div></div>
</div><!-- fragment --><h2><a class="anchor" id="usage_config"></a>
Configuration Management</h2>
<p>The <code>Config</code> singleton manages settings from a YAML configuration file.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="config_8h.html">config.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> Config&amp; config = Config::getInstance();</div>
<div class="line"> <span class="keywordflow">if</span> (config.loadConfig(<span class="stringliteral">&quot;path/to/your/config.yaml&quot;</span>)) {</div>
<div class="line"> <span class="comment">// Get a string value, with a default if not found</span></div>
<div class="line"> std::string outputPath = config.get&lt;std::string&gt;(<span class="stringliteral">&quot;Output:Path&quot;</span>, <span class="stringliteral">&quot;./output/&quot;</span>);</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Output path: &quot;</span> &lt;&lt; outputPath &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get an integer value</span></div>
<div class="line"> <span class="keywordtype">int</span> maxIter = config.get&lt;<span class="keywordtype">int</span>&gt;(<span class="stringliteral">&quot;Solver:MaxIterations&quot;</span>, 100);</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Max iterations: &quot;</span> &lt;&lt; maxIter &lt;&lt; std::endl;</div>
<div class="line"> } <span class="keywordflow">else</span> {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Failed to load configuration.&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
</div><!-- fragment --><h2><a class="anchor" id="usage_logging"></a>
Logging</h2>
<p>The <code>Probe::LogManager</code> provides a way to manage and use loggers.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="probe_8h.html">probe.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="config_8h.html">config.h</a>&quot;</span> <span class="comment">// Often used to configure logging</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> <span class="comment">// Assuming config is loaded and might define log file, level etc.</span></div>
<div class="line"> <span class="comment">// Config::getInstance().loadConfig(&quot;config.yaml&quot;);</span></div>
<div class="line"> </div>
<div class="line"> Probe::LogManager&amp; logManager = Probe::LogManager::getInstance();</div>
<div class="line"> quill::Logger* mainLogger = logManager.getLogger(<span class="stringliteral">&quot;main_app_log&quot;</span>); <span class="comment">// Get or create logger</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Example: Create a new file logger if not configured through a central mechanism</span></div>
<div class="line"> <span class="comment">// quill::Logger* fileLogger = logManager.newFileLogger(&quot;app_trace.log&quot;, &quot;trace_log&quot;);</span></div>
<div class="line"> </div>
<div class="line"> LOG_INFO(mainLogger, <span class="stringliteral">&quot;Application started. Version: {}&quot;</span>, <span class="stringliteral">&quot;1.0.0&quot;</span>);</div>
<div class="line"> <span class="comment">// ... application logic ...</span></div>
<div class="line"> LOG_ERROR(mainLogger, <span class="stringliteral">&quot;An unexpected error occurred in module X.&quot;</span>);</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
</div><!-- fragment --><h2><a class="anchor" id="usage_eos"></a>
Equation of State (EOS)</h2>
<p>The <code>EosIO</code> class loads EOS tables, and the <code>helmholtz</code> namespace provides functions to use them, for example, the Helmholtz EOS.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="_e_o_sio_8h.html">eosIO.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;<a class="code" href="helm_8h.html">helm.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> <span class="keywordflow">try</span> {</div>
<div class="line"> <span class="comment">// Load the Helmholtz EOS table</span></div>
<div class="line"> EosIO helm_eos_io(<span class="stringliteral">&quot;path/to/helm_table.dat&quot;</span>); <span class="comment">// Replace with actual path</span></div>
<div class="line"> EOSTable&amp; table_variant = helm_eos_io.getTable();</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Assuming it&#39;s a HELMTable, get it (add error checking in real code)</span></div>
<div class="line"> <span class="keyword">auto</span>* helm_table_ptr = std::get_if&lt;std::unique_ptr&lt;helmholtz::HELMTable&gt;&gt;(&amp;table_variant);</div>
<div class="line"> <span class="keywordflow">if</span> (!helm_table_ptr || !(*helm_table_ptr) || !(*helm_table_ptr)-&gt;loaded) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Failed to load or access HELM table.&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> <span class="keywordflow">return</span> 1;</div>
<div class="line"> }</div>
<div class="line"> <span class="keyword">const</span> helmholtz::HELMTable&amp; table = **helm_table_ptr;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Define input conditions</span></div>
<div class="line"> helmholtz::EOSInput input;</div>
<div class="line"> input.T = 1.0e7; <span class="comment">// Temperature in K</span></div>
<div class="line"> input.rho = 100.0; <span class="comment">// Density in g/cm^3</span></div>
<div class="line"> input.abar = 1.0; <span class="comment">// Mean atomic mass (e.g., for pure hydrogen)</span></div>
<div class="line"> input.zbar = 1.0; <span class="comment">// Mean atomic number (e.g., for pure hydrogen)</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get EOS results</span></div>
<div class="line"> helmholtz::EOS results = helmholtz::get_helm_EOS(input, table);</div>
<div class="line"> </div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Total Pressure (Ptot): &quot;</span> &lt;&lt; results.ptot &lt;&lt; <span class="stringliteral">&quot; dyne/cm^2&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Total Energy (Etot): &quot;</span> &lt;&lt; results.etot &lt;&lt; <span class="stringliteral">&quot; erg/g&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> } <span class="keywordflow">catch</span> (<span class="keyword">const</span> std::exception&amp; e) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;EOS error: &quot;</span> &lt;&lt; e.what() &lt;&lt; std::endl;</div>
<div class="line"> <span class="keywordflow">return</span> 1;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="a_e_o_sio_8h_html"><div class="ttname"><a href="_e_o_sio_8h.html">EOSio.h</a></div></div>
<div class="ttc" id="ahelm_8h_html"><div class="ttname"><a href="helm_8h.html">helm.h</a></div></div>
</div><!-- fragment --><h2><a class="anchor" id="usage_meshio"></a>
Mesh Handling</h2>
<p>The <code>MeshIO</code> class facilitates loading and managing computational meshes.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &quot;<a class="code" href="mesh_i_o_8h.html">meshIO.h</a>&quot;</span></div>
<div class="line"><span class="preprocessor">#include &quot;mfem.hpp&quot;</span> <span class="comment">// For mfem::Mesh</span></div>
<div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> <a class="code hl_function" href="comp_8cpp.html#ac4c0f8a8146b128f1b8f920e3a9c3b1e">main</a>() {</div>
<div class="line"> <span class="keywordflow">try</span> {</div>
<div class="line"> <span class="comment">// Initialize MeshIO with a mesh file and a scale factor</span></div>
<div class="line"> MeshIO mesh_handler(<span class="stringliteral">&quot;path/to/your/mesh.msh&quot;</span>, 1.0); <span class="comment">// Replace with actual path</span></div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">if</span> (mesh_handler.IsLoaded()) {</div>
<div class="line"> mfem::Mesh&amp; mesh = mesh_handler.GetMesh();</div>
<div class="line"> std::cout &lt;&lt; <span class="stringliteral">&quot;Mesh loaded successfully with &quot;</span> &lt;&lt; mesh.GetNE() &lt;&lt; <span class="stringliteral">&quot; elements.&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Optionally, rescale the mesh</span></div>
<div class="line"> <span class="comment">// mesh_handler.LinearRescale(2.0);</span></div>
<div class="line"> <span class="comment">// std::cout &lt;&lt; &quot;Mesh rescaled. New bounding box: &quot;;</span></div>
<div class="line"> <span class="comment">// mfem::Vector min, max;</span></div>
<div class="line"> <span class="comment">// mesh.GetBoundingBox(min, max);</span></div>
<div class="line"> <span class="comment">// min.Print(std::cout); max.Print(std::cout);</span></div>
<div class="line"> } <span class="keywordflow">else</span> {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;Failed to load mesh.&quot;</span> &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> } <span class="keywordflow">catch</span> (<span class="keyword">const</span> std::exception&amp; e) {</div>
<div class="line"> std::cerr &lt;&lt; <span class="stringliteral">&quot;MeshIO error: &quot;</span> &lt;&lt; e.what() &lt;&lt; std::endl;</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="amesh_i_o_8h_html"><div class="ttname"><a href="mesh_i_o_8h.html">meshIO.h</a></div></div>
</div><!-- fragment --><h1><a class="anchor" id="modules_sec"></a>
Key Modules and Components</h1>
<p>4DSSE is organized into several key modules:</p>
<ul>
<li><b>Polytrope Solver (<code><a class="el" href="poly_solver_8h.html">polySolver.h</a></code>, <code><a class="el" href="polytrope_operator_8h.html">polytropeOperator.h</a></code>):</b> Provides tools to solve the Lane-Emden equation for polytropic stellar structures using a mixed finite element method. <code>PolytropeOperator</code> defines the nonlinear system and its Jacobian, while <code>PolySolver</code> orchestrates the solution process. The <code>SchurCompliment</code> and <code>GMRESInverter</code> classes are helper components for the linear algebra involved.</li>
<li><b>Equation of State (EOS) (<code><a class="el" href="helm_8h.html">helm.h</a></code>, <code><a class="el" href="_e_o_sio_8h.html">eosIO.h</a></code>):</b> Manages Equation of State data. <code><a class="el" href="helm_8h.html">helm.h</a></code> provides an implementation of the Helmholtz EOS (Timmes &amp; Swesty 2000), including structures for table data (<code>HELMTable</code>), input parameters (<code>EOSInput</code>), and output results (<code>EOS</code>). It also defines functions for reading tables and calculating EOS quantities. <code><a class="el" href="_e_o_sio_8h.html">eosIO.h</a></code> provides the <code>EosIO</code> class for loading EOS tables from files, currently supporting the HELM table format.</li>
<li><b>Chemical Composition (<code><a class="el" href="composition_8h.html">composition.h</a></code>, <code>atomicSpecies.h</code>):</b> Manages chemical compositions, allowing representation in mass or number fractions. It interfaces with <code>atomicSpecies.h</code> which provides a database of atomic species properties (based on AME2020).</li>
<li><b>Nuclear Reaction Networks (<code><a class="el" href="network_8h.html">network.h</a></code>, <code><a class="el" href="approx8_8h.html" title="Header file for the Approx8 nuclear reaction network.">approx8.h</a></code>):</b> Defines a base <code>Network</code> class for nuclear reaction network calculations. <code><a class="el" href="approx8_8h.html" title="Header file for the Approx8 nuclear reaction network.">approx8.h</a></code> provides a specific implementation, <code>Approx8Network</code>, for an 8-isotope network (H1, He3, He4, C12, N14, O16, Ne20, Mg24) based on Frank Timmes' "aprox8". It includes functions for individual reaction rates and uses Boost.Numeric.Odeint for solving the ODE system.</li>
<li><b>Physical Constants (<code><a class="el" href="const_8h.html">const.h</a></code>):</b> A singleton class <code>Constants</code> that loads and provides access to a wide range of physical constants with their values, uncertainties, units, and references.</li>
<li><b>Configuration Management (<code><a class="el" href="config_8h.html">config.h</a></code>):</b> A singleton class <code>Config</code> for loading and accessing application settings from YAML configuration files.</li>
<li><b>Probing and Logging (<code><a class="el" href="probe_8h.html">probe.h</a></code>):</b> The <code>Probe</code> namespace offers utility functions for debugging, such as GLVis visualization (<code>glVisView</code>), and a <code>LogManager</code> for handling application-wide logging using the Quill library.</li>
<li><b>Mesh I/O (<code><a class="el" href="mesh_i_o_8h.html">meshIO.h</a></code>):</b> The <code>MeshIO</code> class handles loading and basic manipulation (e.g., scaling) of computational meshes using MFEM's <code>mfem::Mesh</code>. It ensures that meshes are correctly loaded and accessible.</li>
<li><b>Integrators (<code><a class="el" href="integrators_8h.html" title="A collection of utilities for working with MFEM and solving the lane-emden equation.">integrators.h</a></code>):</b> (Details inferred) Likely contains custom MFEM integrators or coefficients used in the finite element formulations.</li>
<li><b>Custom Types (<code><a class="el" href="4_d_s_t_a_r_types_8h.html">4DSTARTypes.h</a></code>):</b> Defines project-specific data type aliases within the <code>SSE</code> namespace, primarily for simplifying common <code>std::pair</code> combinations involving <code>mfem::Array&lt;int&gt;</code> and <code>mfem::Array&lt;double&gt;</code>. These include <code>SSE::MFEMArrayPair</code> and <code>SSE::MFEMArrayPairSet</code>, often used for managing collections of MFEM degree-of-freedom lists and their corresponding values, especially for boundary conditions.</li>
</ul>
<h1><a class="anchor" id="future_dev"></a>
Future Development</h1>
<p>Future work will focus on expanding the physics modules (e.g., equation of state, opacity), improving numerical solvers, and enhancing the parallelization capabilities for large-scale simulations.</p>
<h1><a class="anchor" id="contact_sec"></a>
Contact and Contributions</h1>
<p>For questions, bug reports, or contributions, please refer to the project's repository or contact the development team. (Details to be added) </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>