SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
approx8.h
Go to the documentation of this file.
1/* ***********************************************************************
2//
3// Copyright (C) 2025 -- The 4D-STAR Collaboration
4// File Author: Emily Boudreaux
5// Last Modified: March 21, 2025
6//
7// 4DSSE is free software; you can use it and/or modify
8// it under the terms and restrictions the GNU General Library Public
9// License version 3 (GPLv3) as published by the Free Software Foundation.
10//
11// 4DSSE is distributed in the hope that it will be useful,
12// but WITHOUT ANY WARRANTY; without even the implied warranty of
13// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14// See the GNU Library General Public License for more details.
15//
16// You should have received a copy of the GNU Library General Public License
17// along with this software; if not, write to the Free Software
18// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19//
20// *********************************************************************** */
21#pragma once
22
23#include <array>
24
25#include <boost/numeric/odeint.hpp>
26
27#include "network.h"
28
37
38
40
45 typedef boost::numeric::ublas::vector< double > vector_type;
46
51 typedef boost::numeric::ublas::matrix< double > matrix_type;
52
57 typedef std::array<double,7> vec7;
58
63 struct Approx8Net{
64 static constexpr int ih1=0;
65 static constexpr int ihe3=1;
66 static constexpr int ihe4=2;
67 static constexpr int ic12=3;
68 static constexpr int in14=4;
69 static constexpr int io16=5;
70 static constexpr int ine20=6;
71 static constexpr int img24=7;
72
73 static constexpr int iTemp=img24+1;
74 static constexpr int iDensity =iTemp+1;
75 static constexpr int iEnergy=iDensity+1;
76
77 static constexpr int nIso=img24+1; // number of isotopes
78 static constexpr int nVar=iEnergy+1; // number of variables
79
80 static constexpr std::array<int,nIso> aIon = {
81 1,
82 3,
83 4,
84 12,
85 14,
86 16,
87 20,
88 24
89 };
90
91 static constexpr std::array<double,nIso> mIon = {
92 1.67262164e-24,
93 5.00641157e-24,
94 6.64465545e-24,
95 1.99209977e-23,
96 2.32462686e-23,
97 2.65528858e-23,
98 3.31891077e-23,
99 3.98171594e-23
100 };
101
102 };
103
116 double sum_product( const vec7 &a, const vec7 &b);
117
128 vec7 get_T9_array(const double &T);
129
142 double rate_fit(const vec7 &T9, const vec7 &coef);
143
149 double pp_rate(const vec7 &T9);
150
156 double dp_rate(const vec7 &T9);
157
163 double he3he3_rate(const vec7 &T9);
164
170 double he3he4_rate(const vec7 &T9);
171
177 double triple_alpha_rate(const vec7 &T9);
178
184 double c12p_rate(const vec7 &T9);
185
191 double c12a_rate(const vec7 &T9);
192
198 double n14p_rate(const vec7 &T9);
199
205 double n14a_rate(const vec7 &T9);
206
212 double n15pa_rate(const vec7 &T9);
213
219 double n15pg_rate(const vec7 &T9);
220
226 double n15pg_frac(const vec7 &T9);
227
233 double o16p_rate(const vec7 &T9);
234
240 double o16a_rate(const vec7 &T9);
241
247 double ne20a_rate(const vec7 &T9);
248
254 double c12c12_rate(const vec7 &T9);
255
261 double c12o16_rate(const vec7 &T9);
262
267 struct Jacobian {
274 void operator() ( const vector_type &y, matrix_type &J, double /* t */, vector_type &dfdt ) const;
275 };
276
281 struct ODE {
287 void operator() ( const vector_type &y, vector_type &dydt, double /* t */) const;
288 };
289
294 class Approx8Network final : public Network {
295 public:
297
303 NetOut evaluate(const NetIn &netIn) override;
304
309 void setStiff(bool stiff) override;
310
315 bool isStiff() const override { return m_stiff; }
316 private:
318 double m_tMax = 0;
319 double m_dt0 = 0;
320 bool m_stiff = false;
321
327 static vector_type convert_netIn(const NetIn &netIn);
328 };
329
330
331} // namespace nnApprox8
Network(const NetworkFormat format=NetworkFormat::APPROX8)
Definition network.cpp:32
static vector_type convert_netIn(const NetIn &netIn)
Converts the input parameters to the internal state vector.
Definition approx8.cpp:513
NetOut evaluate(const NetIn &netIn) override
Evaluates the nuclear network.
Definition approx8.cpp:449
void setStiff(bool stiff) override
Sets whether the solver should use a stiff method.
Definition approx8.cpp:509
bool isStiff() const override
Checks if the solver is using a stiff method.
Definition approx8.h:315
double ne20a_rate(const vec7 &T9)
Calculates the rate for the reaction ne20(a,g)mg24.
Definition approx8.cpp:222
double n14p_rate(const vec7 &T9)
Calculates the rate for the reaction n14(p,g)o15 - o15 + p -> c12 + he4.
Definition approx8.cpp:168
double dp_rate(const vec7 &T9)
Calculates the rate for the reaction p + d -> he3.
Definition approx8.cpp:126
double rate_fit(const vec7 &T9, const vec7 &coef)
Definition approx8.cpp:114
double pp_rate(const vec7 &T9)
Calculates the rate for the reaction p + p -> d.
Definition approx8.cpp:119
double o16a_rate(const vec7 &T9)
Calculates the rate for the reaction o16(a,g)ne20.
Definition approx8.cpp:214
double o16p_rate(const vec7 &T9)
Calculates the rate for the reaction o16(p,g)f17 then f17 -> o17(p,a)n14.
Definition approx8.cpp:208
double n15pg_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,g)o16 (CNO II).
Definition approx8.cpp:194
double he3he4_rate(const vec7 &T9)
Calculates the rate for the reaction he3(he3,2p)he4.
Definition approx8.cpp:139
double c12c12_rate(const vec7 &T9)
Calculates the rate for the reaction c12(c12,a)ne20.
Definition approx8.cpp:231
boost::numeric::ublas::matrix< double > matrix_type
Alias for a matrix of doubles using Boost uBLAS.
Definition approx8.h:51
double c12o16_rate(const vec7 &T9)
Calculates the rate for the reaction c12(o16,a)mg24.
Definition approx8.cpp:237
boost::numeric::ublas::vector< double > vector_type
Alias for a vector of doubles using Boost uBLAS.
Definition approx8.h:45
double n15pa_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,a)c12 (CNO I).
Definition approx8.cpp:185
double c12a_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + he4 -> o16.
Definition approx8.cpp:161
vec7 get_T9_array(const double &T)
Definition approx8.cpp:97
double sum_product(const vec7 &a, const vec7 &b)
Definition approx8.cpp:82
double triple_alpha_rate(const vec7 &T9)
Calculates the rate for the reaction he4 + he4 + he4 -> c12.
Definition approx8.cpp:146
double n14a_rate(const vec7 &T9)
Calculates the rate for the reaction n14(a,g)f18 assumed to go on to ne20.
Definition approx8.cpp:177
double he3he3_rate(const vec7 &T9)
Calculates the rate for the reaction he3 + he3 -> he4 + 2p.
Definition approx8.cpp:133
double c12p_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + p -> n13.
Definition approx8.cpp:154
double n15pg_frac(const vec7 &T9)
Calculates the fraction for the reaction n15(p,g)o16.
Definition approx8.cpp:201
std::array< double, 7 > vec7
Alias for a std::array of 7 doubles.
Definition approx8.h:57
Input structure for the network evaluation.
Definition network.h:65
Output structure for the network evaluation.
Definition network.h:89
Contains constants and arrays related to the nuclear network.
Definition approx8.h:63
static constexpr int ic12
Definition approx8.h:67
static constexpr int iTemp
Definition approx8.h:73
static constexpr int ihe4
Definition approx8.h:66
static constexpr int io16
Definition approx8.h:69
static constexpr int nIso
Definition approx8.h:77
static constexpr int ih1
Definition approx8.h:64
static constexpr std::array< int, nIso > aIon
Definition approx8.h:80
static constexpr int img24
Definition approx8.h:71
static constexpr std::array< double, nIso > mIon
Definition approx8.h:91
static constexpr int iEnergy
Definition approx8.h:75
static constexpr int ihe3
Definition approx8.h:65
static constexpr int iDensity
Definition approx8.h:74
static constexpr int in14
Definition approx8.h:68
static constexpr int ine20
Definition approx8.h:70
static constexpr int nVar
Definition approx8.h:78
Functor to calculate the Jacobian matrix for implicit solvers.
Definition approx8.h:267
void operator()(const vector_type &y, matrix_type &J, double, vector_type &dfdt) const
Calculates the Jacobian matrix.
Definition approx8.cpp:247
Functor to calculate the derivatives for the ODE solver.
Definition approx8.h:281
void operator()(const vector_type &y, vector_type &dydt, double) const
Calculates the derivatives of the state vector.
Definition approx8.cpp:352