Multiplies two arrays and sums the resulting elements.
Multiplies two arrays and sums the resulting elements.
#pragma once
#include <array>
#include <boost/numeric/odeint.hpp>
typedef boost::numeric::ublas::vector< double >
vector_type;
typedef boost::numeric::ublas::matrix< double >
matrix_type;
typedef std::array<double,7>
vec7;
struct Approx8Net{
static constexpr int ih1=0;
static constexpr int ihe3=1;
static constexpr int ihe4=2;
static constexpr int ic12=3;
static constexpr int in14=4;
static constexpr int io16=5;
static constexpr int ine20=6;
static constexpr int img24=7;
static constexpr std::array<int,nIso>
aIon = {
1,
3,
4,
12,
14,
16,
20,
24
};
static constexpr std::array<double,nIso>
mIon = {
1.67262164e-24,
5.00641157e-24,
6.64465545e-24,
1.99209977e-23,
2.32462686e-23,
2.65528858e-23,
3.31891077e-23,
3.98171594e-23
};
};
struct Jacobian {
};
struct ODE {
};
class Approx8Network final : public Network {
public:
NetOut
evaluate(
const NetIn &netIn)
override;
private:
};
}
static vector_type convert_netIn(const NetIn &netIn)
Converts the input parameters to the internal state vector.
Definition engine_approx8.cpp:509
bool isStiff() const override
Checks if the solver is using a stiff method.
Definition engine_approx8.h:315
Approx8Network()
Definition engine_approx8.cpp:443
bool m_stiff
Definition engine_approx8.h:320
double m_dt0
Definition engine_approx8.h:319
double m_tMax
Definition engine_approx8.h:318
NetOut evaluate(const NetIn &netIn) override
Evaluates the nuclear network.
Definition engine_approx8.cpp:445
vector_type m_y
Definition engine_approx8.h:317
void setStiff(bool stiff) override
Sets whether the solver should use a stiff method.
Definition engine_approx8.cpp:505
Definition engine_approx8.h:39
double he3he3_rate(const vec7 &T9)
Calculates the rate for the reaction he3 + he3 -> he4 + 2p.
Definition engine_approx8.cpp:129
double pp_rate(const vec7 &T9)
Calculates the rate for the reaction p + p -> d.
Definition engine_approx8.cpp:115
vec7 get_T9_array(const double &T)
Definition engine_approx8.cpp:93
double triple_alpha_rate(const vec7 &T9)
Calculates the rate for the reaction he4 + he4 + he4 -> c12.
Definition engine_approx8.cpp:142
boost::numeric::ublas::matrix< double > matrix_type
Alias for a matrix of doubles using Boost uBLAS.
Definition engine_approx8.h:51
double n14p_rate(const vec7 &T9)
Calculates the rate for the reaction n14(p,g)o15 - o15 + p -> c12 + he4.
Definition engine_approx8.cpp:164
double n14a_rate(const vec7 &T9)
Calculates the rate for the reaction n14(a,g)f18 assumed to go on to ne20.
Definition engine_approx8.cpp:173
double dp_rate(const vec7 &T9)
Calculates the rate for the reaction p + d -> he3.
Definition engine_approx8.cpp:122
double he3he4_rate(const vec7 &T9)
Calculates the rate for the reaction he3(he3,2p)he4.
Definition engine_approx8.cpp:135
double o16p_rate(const vec7 &T9)
Calculates the rate for the reaction o16(p,g)f17 then f17 -> o17(p,a)n14.
Definition engine_approx8.cpp:204
double c12c12_rate(const vec7 &T9)
Calculates the rate for the reaction c12(c12,a)ne20.
Definition engine_approx8.cpp:227
double o16a_rate(const vec7 &T9)
Calculates the rate for the reaction o16(a,g)ne20.
Definition engine_approx8.cpp:210
double c12p_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + p -> n13.
Definition engine_approx8.cpp:150
double c12o16_rate(const vec7 &T9)
Calculates the rate for the reaction c12(o16,a)mg24.
Definition engine_approx8.cpp:233
double n15pa_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,a)c12 (CNO I).
Definition engine_approx8.cpp:181
boost::numeric::ublas::vector< double > vector_type
Alias for a vector of doubles using Boost uBLAS.
Definition engine_approx8.h:45
std::array< double, 7 > vec7
Alias for a std::array of 7 doubles.
Definition engine_approx8.h:57
double n15pg_frac(const vec7 &T9)
Calculates the fraction for the reaction n15(p,g)o16.
Definition engine_approx8.cpp:197
double n15pg_rate(const vec7 &T9)
Calculates the rate for the reaction n15(p,g)o16 (CNO II).
Definition engine_approx8.cpp:190
double ne20a_rate(const vec7 &T9)
Calculates the rate for the reaction ne20(a,g)mg24.
Definition engine_approx8.cpp:218
double rate_fit(const vec7 &T9, const vec7 &coef)
Definition engine_approx8.cpp:110
double c12a_rate(const vec7 &T9)
Calculates the rate for the reaction c12 + he4 -> o16.
Definition engine_approx8.cpp:157
static constexpr int iTemp
Definition engine_approx8.h:73
static constexpr int iEnergy
Definition engine_approx8.h:75
static constexpr int in14
Definition engine_approx8.h:68
static constexpr std::array< int, nIso > aIon
Definition engine_approx8.h:80
static constexpr int nIso
Definition engine_approx8.h:77
static constexpr int iDensity
Definition engine_approx8.h:74
static constexpr int nVar
Definition engine_approx8.h:78
static constexpr int ihe4
Definition engine_approx8.h:66
static constexpr std::array< double, nIso > mIon
Definition engine_approx8.h:91
static constexpr int ic12
Definition engine_approx8.h:67
static constexpr int img24
Definition engine_approx8.h:71
static constexpr int ihe3
Definition engine_approx8.h:65
static constexpr int io16
Definition engine_approx8.h:69
static constexpr int ih1
Definition engine_approx8.h:64
static constexpr int ine20
Definition engine_approx8.h:70
void operator()(const vector_type &y, matrix_type &J, double, vector_type &dfdt) const
Calculates the Jacobian matrix.
Definition engine_approx8.cpp:243
void operator()(const vector_type &y, vector_type &dydt, double) const
Calculates the derivatives of the state vector.
Definition engine_approx8.cpp:348