25#include <boost/numeric/odeint.hpp>
27#include "fourdst/constants/const.h"
28#include "fourdst/config/config.h"
29#include "quill/LogMacros.h"
78 using namespace boost::numeric::odeint;
84 for (
size_t i=0; i < a.size(); i++) {
95 const double T9=1e-9*T;
96 const double T913=pow(T9,1./3.);
103 arr[5]=pow(T9,5./3.);
116 constexpr vec7 a1 = {-34.78630, 0,-3.511930, 3.100860, -0.1983140, 1.262510e-2, -1.025170};
117 constexpr vec7 a2 = { -4.364990e+1,-2.460640e-3,-2.750700,-4.248770e-1,1.598700e-2,-6.908750e-4,-2.076250e-1};
123 constexpr vec7 a1 = {7.528980, 0, -3.720800, 0.8717820, 0, 0,-0.6666670};
124 constexpr vec7 a2 = {8.935250, 0, -3.720800, 0.1986540, 0, 0, 0.3333330};
130 constexpr vec7 a = {2.477880e+01,0,-12.27700,-0.1036990,-6.499670e-02,1.681910e-02,-6.666670e-01};
136 constexpr vec7 a1 = {1.560990e+01,0.000000e+00,-1.282710e+01,-3.082250e-02,-6.546850e-01,8.963310e-02,-6.666670e-01};
137 constexpr vec7 a2 = {1.770750e+01,0.000000e+00,-1.282710e+01,-3.812600e+00,9.422850e-02,-3.010180e-03,1.333330e+00};
143 constexpr vec7 a1 = {-9.710520e-01,0.000000e+00,-3.706000e+01,2.934930e+01,-1.155070e+02,-1.000000e+01,-1.333330e+00};
144 constexpr vec7 a2 = {-1.178840e+01,-1.024460e+00,-2.357000e+01,2.048860e+01,-1.298820e+01,-2.000000e+01,-2.166670e+00};
145 constexpr vec7 a3 = {-2.435050e+01,-4.126560e+00,-1.349000e+01,2.142590e+01,-1.347690e+00,8.798160e-02,-1.316530e+01};
151 constexpr vec7 a1={1.714820e+01,0.000000e+00,-1.369200e+01,-2.308810e-01,4.443620e+00,-3.158980e+00,-6.666670e-01};
152 constexpr vec7 a2={1.754280e+01,-3.778490e+00,-5.107350e+00,-2.241110e+00,1.488830e-01,0.000000e+00,-1.500000e+00};
158 constexpr vec7 a1={6.965260e+01,-1.392540e+00,5.891280e+01,-1.482730e+02,9.083240e+00,-5.410410e-01,7.035540e+01};
159 constexpr vec7 a2={2.546340e+02,-1.840970e+00,1.034110e+02,-4.205670e+02,6.408740e+01,-1.246240e+01,1.373030e+02};
165 constexpr vec7 a1={1.701000e+01,0.000000e+00,-1.519300e+01,-1.619540e-01,-7.521230e+00,-9.875650e-01,-6.666670e-01};
166 constexpr vec7 a2={2.011690e+01,0.000000e+00,-1.519300e+01,-4.639750e+00,9.734580e+00,-9.550510e+00,3.333330e-01};
167 constexpr vec7 a3={7.654440e+00,-2.998000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
168 constexpr vec7 a4={6.735780e+00,-4.891000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,6.820000e-02};
174 constexpr vec7 a1={2.153390e+01,0.000000e+00,-3.625040e+01,0.000000e+00,0.000000e+00,-5.000000e+00,-6.666670e-01};
175 constexpr vec7 a2={1.968380e-01,-5.160340e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
176 constexpr vec7 a3={1.389950e+01,-1.096560e+01,-5.622700e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
182 constexpr vec7 a1 = {2.747640e+01,0.000000e+00,-1.525300e+01,1.593180e+00,2.447900e+00,-2.197080e+00,-6.666670e-01};
183 constexpr vec7 a2 = {-4.873470e+00,-2.021170e+00,0.000000e+00,3.084970e+01,-8.504330e+00,-1.544260e+00,-1.500000e+00};
184 constexpr vec7 a3 = {2.089720e+01,-7.406000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
185 constexpr vec7 a4 = {-6.575220e+00,-1.163800e+00,0.000000e+00,2.271050e+01,-2.907070e+00,2.057540e-01,-1.500000e+00};
191 constexpr vec7 a1 = {2.001760e+01,0.000000e+00,-1.524000e+01,3.349260e-01,4.590880e+00,-4.784680e+00,-6.666670e-01};
192 constexpr vec7 a2 = {6.590560e+00,-2.923150e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
193 constexpr vec7 a3 = {1.454440e+01,-1.022950e+01,0.000000e+00,0.000000e+00,4.590370e-02,0.000000e+00,-1.500000e+00};
205 constexpr vec7 a={1.909040e+01,0.000000e+00,-1.669600e+01,-1.162520e+00,2.677030e-01,-3.384110e-02,-6.666670e-01};
211 constexpr vec7 a1={2.390300e+01,0.000000e+00,-3.972620e+01,-2.107990e-01,4.428790e-01,-7.977530e-02,-6.666670e-01};
212 constexpr vec7 a2={3.885710e+00,-1.035850e+01,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
213 constexpr vec7 a3={9.508480e+00,-1.276430e+01,0.000000e+00,-3.659250e+00,7.142240e-01,-1.075080e-03,-1.500000e+00};
219 constexpr vec7 a1={2.450580e+01,0.000000e+00,-4.625250e+01,5.589010e+00,7.618430e+00,-3.683000e+00,-6.666670e-01};
220 constexpr vec7 a2={-3.870550e+01,-2.506050e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
221 constexpr vec7 a3={1.983070e+00,-9.220260e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
222 constexpr vec7 a4={-8.798270e+00,-1.278090e+01,0.000000e+00,1.692290e+01,-2.573250e+00,2.089970e-01,-1.500000e+00};
228 constexpr vec7 a={6.128630e+01,0.000000e+00,-8.416500e+01,-1.566270e+00,-7.360840e-02,-7.279700e-02,-6.666670e-01};
234 constexpr vec7 a={4.853410e+01,3.720400e-01,-1.334130e+02,5.015720e+01,-3.159870e+00,1.782510e-02,-2.370270e+01};
244 fourdst::constant::Constants& constants = fourdst::constant::Constants::getInstance();
245 const double avo = constants.get(
"N_a").value;
246 const double clight = constants.get(
"c").value;
266 const double aFrac=1-pFrac;
349 const fourdst::constant::Constants& constants = fourdst::constant::Constants::getInstance();
350 const double avo = constants.get(
"N_a").value;
351 const double clight = constants.get(
"c").value;
359 const double rpp=den*
pp_rate(T9);
374 const double aFrac=1-pFrac;
376 const double yh1 = y[Approx8Net:: ih1];
450 const double stiff_abs_tol =
m_config.get<
double>(
"Network:Approx8:Stiff:AbsTol", 1.0e-6);
451 const double stiff_rel_tol =
m_config.get<
double>(
"Network:Approx8:Stiff:RelTol", 1.0e-6);
452 const double nonstiff_abs_tol =
m_config.get<
double>(
"Network:Approx8:NonStiff:AbsTol", 1.0e-6);
453 const double nonstiff_rel_tol =
m_config.get<
double>(
"Network:Approx8:NonStiff:RelTol", 1.0e-6);
458 LOG_DEBUG(
m_logger,
"Using stiff solver for Approx8Network");
459 num_steps = integrate_const(
460 make_dense_output<rosenbrock4<double>>(stiff_abs_tol, stiff_rel_tol),
469 LOG_DEBUG(
m_logger,
"Using non stiff solver for Approx8Network");
470 num_steps = integrate_const (
471 make_dense_output<runge_kutta_dopri5<vector_type>>(nonstiff_abs_tol, nonstiff_rel_tol),
490 std::vector<double> outComposition;
494 outComposition.push_back(
m_y[i]);
499 const std::vector<std::string> symbols = {
"H-1",
"He-3",
"He-4",
"C-12",
"N-14",
"O-16",
"Ne-20",
"Mg-24"};
500 netOut.
composition = fourdst::composition::Composition(symbols, outComposition);
Network(const NetworkFormat format=NetworkFormat::APPROX8)
quill::Logger * m_logger
Logger instance.
fourdst::config::Config & m_config
Configuration instance.
static vector_type convert_netIn(const NetIn &netIn)
Converts the input parameters to the internal state vector.
NetOut evaluate(const NetIn &netIn) override
Evaluates the nuclear network.
void setStiff(bool stiff) override
Sets whether the solver should use a stiff method.
double he3he3_rate(const vec7 &T9)
Calculates the rate for the reaction he3 + he3 -> he4 + 2p.
double pp_rate(const vec7 &T9)
Calculates the rate for the reaction p + p -> d.
vec7 get_T9_array(const double &T)
double triple_alpha_rate(const vec7 &T9)
Calculates the rate for the reaction he4 + he4 + he4 -> c12.
boost::numeric::ublas::matrix< double > matrix_type
Alias for a matrix of doubles using Boost uBLAS.
double n14p_rate(const vec7 &T9)
Calculates the rate for the reaction n14(p,g)o15 - o15 + p -> c12 + he4.
double n14a_rate(const vec7 &T9)
Calculates the rate for the reaction n14(a,g)f18 assumed to go on to ne20.
double dp_rate(const vec7 &T9)
Calculates the rate for the reaction p + d -> he3.
double he3he4_rate(const vec7 &T9)
Calculates the rate for the reaction he3(he3,2p)he4.
double o16p_rate(const vec7 &T9)
Calculates the rate for the reaction o16(p,g)f17 then f17 -> o17(p,a)n14.
double c12c12_rate(const vec7 &T9)
Calculates the rate for the reaction c12(c12,a)ne20.
double o16a_rate(const vec7 &T9)
Calculates the rate for the reaction o16(a,g)ne20.
double c12p_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + p -> n13.
double c12o16_rate(const vec7 &T9)
Calculates the rate for the reaction c12(o16,a)mg24.
double n15pa_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,a)c12 (CNO I).
boost::numeric::ublas::vector< double > vector_type
Alias for a vector of doubles using Boost uBLAS.
std::array< double, 7 > vec7
Alias for a std::array of 7 doubles.
double sum_product(const vec7 &a, const vec7 &b)
double n15pg_frac(const vec7 &T9)
Calculates the fraction for the reaction n15(p,g)o16.
double n15pg_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,g)o16 (CNO II).
double ne20a_rate(const vec7 &T9)
Calculates the rate for the reaction ne20(a,g)mg24.
double rate_fit(const vec7 &T9, const vec7 &coef)
double c12a_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + he4 -> o16.
@ APPROX8
Approx8 nuclear reaction network format.
double density
Density in g/cm^3.
fourdst::composition::Composition composition
Composition of the network.
double dt0
Initial time step.
double temperature
Temperature in Kelvin.
double energy
Energy in ergs.
fourdst::composition::Composition composition
Composition of the network after evaluation.
double energy
Energy in ergs after evaluation.
int num_steps
Number of steps taken in the evaluation.
static constexpr int iTemp
static constexpr int iEnergy
static constexpr int in14
static constexpr std::array< int, nIso > aIon
static constexpr int nIso
static constexpr int iDensity
static constexpr int nVar
static constexpr int ihe4
static constexpr std::array< double, nIso > mIon
static constexpr int ic12
static constexpr int img24
static constexpr int ihe3
static constexpr int io16
static constexpr int ine20
Functor to calculate the Jacobian matrix for implicit solvers.
void operator()(const vector_type &y, matrix_type &J, double, vector_type &dfdt) const
Calculates the Jacobian matrix.
Functor to calculate the derivatives for the ODE solver.
void operator()(const vector_type &y, vector_type &dydt, double) const
Calculates the derivatives of the state vector.