SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
approx8.cpp
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#include <cmath>
22#include <stdexcept>
23#include <array>
24
25#include <boost/numeric/odeint.hpp>
26
27#include "const.h"
28#include "config.h"
29#include "quill/LogMacros.h"
30
31#include "approx8.h"
32#include "network.h"
33
34/* Nuclear reaction network in cgs units based on Frank Timmes' "approx8".
35 At this time it does neither screening nor neutrino losses. It includes
36 the following 8 isotopes:
37
38 h1
39 he3
40 he4
41 c12
42 n14
43 o16
44 ne20
45 mg24
46
47 and the following nuclear reactions:
48
49 ---pp chain---
50 p(p,e+)d
51 d(p,g)he3
52 he3(he3,2p)he4
53
54 ---CNO cycle---
55 c12(p,g)n13 - n13 -> c13 + p -> n14
56 n14(p,g)o15 - o15 + p -> c12 + he4
57 n14(a,g)f18 - proceeds to ne20
58 n15(p,a)c12 - / these two n15 reactions are \ CNO I
59 n15(p,g)o16 - \ used to calculate a fraction / CNO II
60 o16(p,g)f17 - f17 + e -> o17 and then o17 + p -> n14 + he4
61
62 ---alpha captures---
63 c12(a,g)o16
64 triple alpha
65 o16(a,g)ne20
66 ne20(a,g)mg24
67 c12(c12,a)ne20
68 c12(o16,a)mg24
69
70At present the rates are all evaluated using a fitting function.
71The coefficients to the fit are from reaclib.jinaweb.org .
72
73*/
74
76
77 // using namespace std;
78 using namespace boost::numeric::odeint;
79
80 //helper functions
81 // a function to multiply two arrays and then sum the resulting elements: sum(a*b)
82 double sum_product( const vec7 &a, const vec7 &b){
83 if (a.size() != b.size()) {
84 throw std::runtime_error("Error: array size mismatch in sum_product");
85 }
86
87 double sum=0;
88 for (size_t i=0; i < a.size(); i++) {
89 sum += a[i] * b[i];
90 }
91 return sum;
92 }
93
94 // the fit to the nuclear reaction rates is of the form:
95 // exp( a0 + a1/T9 + a2/T9^(1/3) + a3*T9^(1/3) + a4*T9 + a5*T9^(5/3) + log(T9) )
96 // this function returns an array of the T9 terms in that order, where T9 is the temperatures in GigaKelvin
97 vec7 get_T9_array(const double &T) {
98 vec7 arr;
99 const double T9=1e-9*T;
100 const double T913=pow(T9,1./3.);
101
102 arr[0]=1;
103 arr[1]=1/T9;
104 arr[2]=1/T913;
105 arr[3]=T913;
106 arr[4]=T9;
107 arr[5]=pow(T9,5./3.);
108 arr[6]=log(T9);
109
110 return arr;
111 }
112
113 // this function uses the two preceding functions to evaluate the rate given the T9 array and the coefficients
114 double rate_fit(const vec7 &T9, const vec7 &coef){
115 return exp(sum_product(T9,coef));
116 }
117
118 // p + p -> d; this, like some of the other rates, this is a composite of multiple fits
119 double pp_rate(const vec7 &T9) {
120 constexpr vec7 a1 = {-34.78630, 0,-3.511930, 3.100860, -0.1983140, 1.262510e-2, -1.025170};
121 constexpr vec7 a2 = { -4.364990e+1,-2.460640e-3,-2.750700,-4.248770e-1,1.598700e-2,-6.908750e-4,-2.076250e-1};
122 return rate_fit(T9,a1) + rate_fit(T9,a2);
123 }
124
125 // p + d -> he3
126 double dp_rate(const vec7 &T9) {
127 constexpr vec7 a1 = {7.528980, 0, -3.720800, 0.8717820, 0, 0,-0.6666670};
128 constexpr vec7 a2 = {8.935250, 0, -3.720800, 0.1986540, 0, 0, 0.3333330};
129 return rate_fit(T9,a1) + rate_fit(T9,a2);
130 }
131
132 // he3 + he3 -> he4 + 2p
133 double he3he3_rate(const vec7 &T9){
134 constexpr vec7 a = {2.477880e+01,0,-12.27700,-0.1036990,-6.499670e-02,1.681910e-02,-6.666670e-01};
135 return rate_fit(T9,a);
136 }
137
138 // he3(he3,2p)he4
139 double he3he4_rate(const vec7 &T9){
140 constexpr vec7 a1 = {1.560990e+01,0.000000e+00,-1.282710e+01,-3.082250e-02,-6.546850e-01,8.963310e-02,-6.666670e-01};
141 constexpr vec7 a2 = {1.770750e+01,0.000000e+00,-1.282710e+01,-3.812600e+00,9.422850e-02,-3.010180e-03,1.333330e+00};
142 return rate_fit(T9,a1) + rate_fit(T9,a2);
143 }
144
145 // he4 + he4 + he4 -> c12
146 double triple_alpha_rate(const vec7 &T9){
147 constexpr vec7 a1 = {-9.710520e-01,0.000000e+00,-3.706000e+01,2.934930e+01,-1.155070e+02,-1.000000e+01,-1.333330e+00};
148 constexpr vec7 a2 = {-1.178840e+01,-1.024460e+00,-2.357000e+01,2.048860e+01,-1.298820e+01,-2.000000e+01,-2.166670e+00};
149 constexpr vec7 a3 = {-2.435050e+01,-4.126560e+00,-1.349000e+01,2.142590e+01,-1.347690e+00,8.798160e-02,-1.316530e+01};
150 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
151 }
152
153 // c12 + p -> n13
154 double c12p_rate(const vec7 &T9){
155 constexpr vec7 a1={1.714820e+01,0.000000e+00,-1.369200e+01,-2.308810e-01,4.443620e+00,-3.158980e+00,-6.666670e-01};
156 constexpr vec7 a2={1.754280e+01,-3.778490e+00,-5.107350e+00,-2.241110e+00,1.488830e-01,0.000000e+00,-1.500000e+00};
157 return rate_fit(T9,a1) + rate_fit(T9,a2);
158 }
159
160 // c12 + he4 -> o16
161 double c12a_rate(const vec7 &T9){
162 constexpr vec7 a1={6.965260e+01,-1.392540e+00,5.891280e+01,-1.482730e+02,9.083240e+00,-5.410410e-01,7.035540e+01};
163 constexpr vec7 a2={2.546340e+02,-1.840970e+00,1.034110e+02,-4.205670e+02,6.408740e+01,-1.246240e+01,1.373030e+02};
164 return rate_fit(T9,a1) + rate_fit(T9,a2);
165 }
166
167 // n14(p,g)o15 - o15 + p -> c12 + he4
168 double n14p_rate(const vec7 &T9){
169 constexpr vec7 a1={1.701000e+01,0.000000e+00,-1.519300e+01,-1.619540e-01,-7.521230e+00,-9.875650e-01,-6.666670e-01};
170 constexpr vec7 a2={2.011690e+01,0.000000e+00,-1.519300e+01,-4.639750e+00,9.734580e+00,-9.550510e+00,3.333330e-01};
171 constexpr vec7 a3={7.654440e+00,-2.998000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
172 constexpr vec7 a4={6.735780e+00,-4.891000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,6.820000e-02};
173 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
174 }
175
176 // n14(a,g)f18 assumed to go on to ne20
177 double n14a_rate(const vec7 &T9){
178 constexpr vec7 a1={2.153390e+01,0.000000e+00,-3.625040e+01,0.000000e+00,0.000000e+00,-5.000000e+00,-6.666670e-01};
179 constexpr vec7 a2={1.968380e-01,-5.160340e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
180 constexpr vec7 a3={1.389950e+01,-1.096560e+01,-5.622700e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
181 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
182 }
183
184 // n15(p,a)c12 (CNO I)
185 double n15pa_rate(const vec7 &T9){
186 constexpr vec7 a1 = {2.747640e+01,0.000000e+00,-1.525300e+01,1.593180e+00,2.447900e+00,-2.197080e+00,-6.666670e-01};
187 constexpr vec7 a2 = {-4.873470e+00,-2.021170e+00,0.000000e+00,3.084970e+01,-8.504330e+00,-1.544260e+00,-1.500000e+00};
188 constexpr vec7 a3 = {2.089720e+01,-7.406000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
189 constexpr vec7 a4 = {-6.575220e+00,-1.163800e+00,0.000000e+00,2.271050e+01,-2.907070e+00,2.057540e-01,-1.500000e+00};
190 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
191 }
192
193 // n15(p,g)o16 (CNO II)
194 double n15pg_rate(const vec7 &T9){
195 constexpr vec7 a1 = {2.001760e+01,0.000000e+00,-1.524000e+01,3.349260e-01,4.590880e+00,-4.784680e+00,-6.666670e-01};
196 constexpr vec7 a2 = {6.590560e+00,-2.923150e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
197 constexpr vec7 a3 = {1.454440e+01,-1.022950e+01,0.000000e+00,0.000000e+00,4.590370e-02,0.000000e+00,-1.500000e+00};
198 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
199 }
200
201 double n15pg_frac(const vec7 &T9){
202 const double f1=n15pg_rate(T9);
203 const double f2=n15pa_rate(T9);
204 return f1/(f1+f2);
205 }
206
207 // o16(p,g)f17 then f17 -> o17(p,a)n14
208 double o16p_rate(const vec7 &T9){
209 constexpr vec7 a={1.909040e+01,0.000000e+00,-1.669600e+01,-1.162520e+00,2.677030e-01,-3.384110e-02,-6.666670e-01};
210 return rate_fit(T9,a);
211 }
212
213 // o16(a,g)ne20
214 double o16a_rate(const vec7 &T9){
215 constexpr vec7 a1={2.390300e+01,0.000000e+00,-3.972620e+01,-2.107990e-01,4.428790e-01,-7.977530e-02,-6.666670e-01};
216 constexpr vec7 a2={3.885710e+00,-1.035850e+01,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
217 constexpr vec7 a3={9.508480e+00,-1.276430e+01,0.000000e+00,-3.659250e+00,7.142240e-01,-1.075080e-03,-1.500000e+00};
218 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
219 }
220
221 // ne20(a,g)mg24
222 double ne20a_rate(const vec7 &T9){
223 constexpr vec7 a1={2.450580e+01,0.000000e+00,-4.625250e+01,5.589010e+00,7.618430e+00,-3.683000e+00,-6.666670e-01};
224 constexpr vec7 a2={-3.870550e+01,-2.506050e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
225 constexpr vec7 a3={1.983070e+00,-9.220260e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,-1.500000e+00};
226 constexpr vec7 a4={-8.798270e+00,-1.278090e+01,0.000000e+00,1.692290e+01,-2.573250e+00,2.089970e-01,-1.500000e+00};
227 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
228 }
229
230 // c12(c12,a)ne20
231 double c12c12_rate(const vec7 &T9){
232 constexpr vec7 a={6.128630e+01,0.000000e+00,-8.416500e+01,-1.566270e+00,-7.360840e-02,-7.279700e-02,-6.666670e-01};
233 return rate_fit(T9,a);
234 }
235
236 // c12(o16,a)mg24
237 double c12o16_rate(const vec7 &T9){
238 constexpr vec7 a={4.853410e+01,3.720400e-01,-1.334130e+02,5.015720e+01,-3.159870e+00,1.782510e-02,-2.370270e+01};
239 return rate_fit(T9,a);
240 }
241
242
243 // for Boost ODE solvers either a struct or a class is required
244
245 // a Jacobian matrix for implicit solvers
246
247 void Jacobian::operator() ( const vector_type &y, matrix_type &J, double /* t */, vector_type &dfdt ) const {
249 const double avo = constants.get("N_a").value;
250 const double clight = constants.get("c").value;
251 // EOS
253
254 // evaluate rates once per call
255 const double rpp=pp_rate(T9);
256 const double r33=he3he3_rate(T9);
257 const double r34=he3he4_rate(T9);
258 const double r3a=triple_alpha_rate(T9);
259 const double rc12p=c12p_rate(T9);
260 const double rc12a=c12a_rate(T9);
261 const double rn14p=n14p_rate(T9);
262 const double rn14a=n14a_rate(T9);
263 const double ro16p=o16p_rate(T9);
264 const double ro16a=o16a_rate(T9);
265 const double rne20a=ne20a_rate(T9);
266 const double r1212=c12c12_rate(T9);
267 const double r1216=c12o16_rate(T9);
268
269 const double pFrac=n15pg_frac(T9);
270 const double aFrac=1-pFrac;
271
272 const double yh1 = y[Approx8Net::ih1];
273 const double yhe3 = y[Approx8Net::ihe3];
274 const double yhe4 = y[Approx8Net::ihe4];
275 const double yc12 = y[Approx8Net::ic12];
276 const double yn14 = y[Approx8Net::in14];
277 const double yo16 = y[Approx8Net::io16];
278 const double yne20 = y[Approx8Net::ine20];
279
280 // zero all elements to begin
281 for (int i=0; i < Approx8Net::nVar; i++) {
282 for (int j=0; j < Approx8Net::nVar; j++) {
283 J(i,j)=0.0;
284 }
285 }
286
287 // h1 jacobian elements
288 J(Approx8Net::ih1,Approx8Net::ih1) = -3*yh1*rpp - 2*yc12*rc12p -2*yn14*rn14p -2*yo16*ro16p;
289 J(Approx8Net::ih1,Approx8Net::ihe3) = 2*yhe3*r33 - yhe4*r34;
290 J(Approx8Net::ih1,Approx8Net::ihe4) = -yhe3*r34;
291 J(Approx8Net::ih1,Approx8Net::ic12) = -2*yh1*rc12p;
292 J(Approx8Net::ih1,Approx8Net::in14) = -2*yh1*rn14p;
293 J(Approx8Net::ih1,Approx8Net::io16) = -2*yh1*ro16p;
294
295 // he3 jacobian elements
297 J(Approx8Net::ihe3,Approx8Net::ihe3) = -2*yhe3*r33 - yhe4*r34;
298 J(Approx8Net::ihe3,Approx8Net::ihe4) = -yhe3*r34;
299
300 // he4 jacobian elements
301 J(Approx8Net::ihe4,Approx8Net::ih1) = yn14*aFrac*rn14p + yo16*ro16p;
302 J(Approx8Net::ihe4,Approx8Net::ihe3) = yhe3*r33 - yhe4*r34;
303 J(Approx8Net::ihe4,Approx8Net::ihe4) = yhe3*r34 - 1.5*yhe4*yhe4*r3a - yc12*rc12a - 1.5*yn14*rn14a - yo16*ro16a - yne20*rne20a;
304 J(Approx8Net::ihe4,Approx8Net::ic12) = -yhe4*rc12a + yc12*r1212 + yo16*r1216;
305 J(Approx8Net::ihe4,Approx8Net::in14) = yh1*aFrac*rn14p - 1.5*yhe4*rn14a;
306 J(Approx8Net::ihe4,Approx8Net::io16) = yh1*ro16p - yhe4*ro16a + yc12*r1216;
307 J(Approx8Net::ihe4,Approx8Net::ine20) = -yhe4*rne20a;
308
309 // c12 jacobian elements
310 J(Approx8Net::ic12,Approx8Net::ih1) = -yc12*rc12p + yn14*aFrac*rn14p;
311 J(Approx8Net::ic12,Approx8Net::ihe4) = 0.5*yhe4*yhe4*r3a - yhe4*rc12a;
312 J(Approx8Net::ic12,Approx8Net::ic12) = -yh1*rc12p - yhe4*rc12a - yo16*r1216 - 2*yc12*r1212;
313 J(Approx8Net::ic12,Approx8Net::in14) = yh1*yn14*aFrac*rn14p;
314 J(Approx8Net::ic12,Approx8Net::io16) = -yc12*r1216;
315
316 // n14 jacobian elements
317 J(Approx8Net::in14,Approx8Net::ih1) = yc12*rc12p - yn14*rn14p + yo16*ro16p;
318 J(Approx8Net::in14,Approx8Net::ihe4) = -yn14*rn14a;
319 J(Approx8Net::in14,Approx8Net::ic12) = yh1*rc12p;
320 J(Approx8Net::in14,Approx8Net::in14) = -yh1*rn14p - yhe4*rn14a;
321 J(Approx8Net::in14,Approx8Net::io16) = yo16*ro16p;
322
323 // o16 jacobian elements
324 J(Approx8Net::io16,Approx8Net::ih1) = yn14*pFrac*rn14p - yo16*ro16p;
325 J(Approx8Net::io16,Approx8Net::ihe4) = yc12*rc12a - yo16*ro16a;
326 J(Approx8Net::io16,Approx8Net::ic12) = yhe4*rc12a - yo16*r1216;
327 J(Approx8Net::io16,Approx8Net::in14) = yh1*pFrac*rn14p;
328 J(Approx8Net::io16,Approx8Net::io16) = yh1*ro16p - yc12*r1216 -yhe4*ro16a;
329
330 // ne20 jacobian elements
331 J(Approx8Net::ine20,Approx8Net::ihe4) = yn14*rn14a + yo16*ro16a - yne20*rne20a;
332 J(Approx8Net::ine20,Approx8Net::ic12) = yc12*r1212;
333 J(Approx8Net::ine20,Approx8Net::in14) = yhe4*rn14a;
334 J(Approx8Net::ine20,Approx8Net::io16) = yo16*ro16a;
335 J(Approx8Net::ine20,Approx8Net::ine20) = -yhe4*rne20a;
336
337 // mg24 jacobian elements
338 J(Approx8Net::img24,Approx8Net::ihe4) = yne20*rne20a;
339 J(Approx8Net::img24,Approx8Net::ic12) = yo16*r1216;
340 J(Approx8Net::img24,Approx8Net::io16) = yc12*r1216;
341 J(Approx8Net::img24,Approx8Net::ine20) = yhe4*rne20a;
342
343 // energy accounting
344 for (int j=0; j<Approx8Net::nIso; j++) {
345 for (int i=0; i<Approx8Net::nIso; i++) {
346 J(Approx8Net::iEnergy,j) += J(i,j)*Approx8Net::mIon[i];
347 }
348 J(Approx8Net::iEnergy,j) *= -avo*clight*clight;
349 }
350 }
351
352 void ODE::operator() ( const vector_type &y, vector_type &dydt, double /* t */) const {
354 const double avo = constants.get("N_a").value;
355 const double clight = constants.get("c").value;
356
357 // EOS
358 const double T = y[Approx8Net::iTemp];
359 const double den = y[Approx8Net::iDensity];
360 const vec7 T9=get_T9_array(T);
361
362 // rates
363 const double rpp=den*pp_rate(T9);
364 const double r33=den*he3he3_rate(T9);
365 const double r34=den*he3he4_rate(T9);
366 const double r3a=den*den*triple_alpha_rate(T9);
367 const double rc12p=den*c12p_rate(T9);
368 const double rc12a=den*c12a_rate(T9);
369 const double rn14p=den*n14p_rate(T9);
370 const double rn14a=n14a_rate(T9);
371 const double ro16p=den*o16p_rate(T9);
372 const double ro16a=den*o16a_rate(T9);
373 const double rne20a=den*ne20a_rate(T9);
374 const double r1212=den*c12c12_rate(T9);
375 const double r1216=den*c12o16_rate(T9);
376
377 const double pFrac=n15pg_frac(T9);
378 const double aFrac=1-pFrac;
379
380 const double yh1 = y[Approx8Net:: ih1];
381 const double yhe3 = y[Approx8Net:: ihe3];
382 const double yhe4 = y[Approx8Net:: ihe4];
383 const double yc12 = y[Approx8Net:: ic12];
384 const double yn14 = y[Approx8Net:: in14];
385 const double yo16 = y[Approx8Net:: io16];
386 const double yne20 = y[Approx8Net::ine20];
387
388 dydt[Approx8Net::ih1] = -1.5*yh1*yh1*rpp;
389 dydt[Approx8Net::ih1] += yhe3*yhe3*r33;
390 dydt[Approx8Net::ih1] += -yhe3*yhe4*r34;
391 dydt[Approx8Net::ih1] += -2*yh1*yc12*rc12p;
392 dydt[Approx8Net::ih1] += -2*yh1*yn14*rn14p;
393 dydt[Approx8Net::ih1] += -2*yh1*yo16*ro16p;
394
395 dydt[Approx8Net::ihe3] = 0.5*yh1*yh1*rpp;
396 dydt[Approx8Net::ihe3] += -yhe3*yhe3*r33;
397 dydt[Approx8Net::ihe3] += -yhe3*yhe4*r34;
398
399 dydt[Approx8Net::ihe4] = 0.5*yhe3*yhe3*r33;
400 dydt[Approx8Net::ihe4] += yhe3*yhe4*r34;
401 dydt[Approx8Net::ihe4] += -yhe4*yc12*rc12a;
402 dydt[Approx8Net::ihe4] += yh1*yn14*aFrac*rn14p;
403 dydt[Approx8Net::ihe4] += yh1*yo16*ro16p;
404 dydt[Approx8Net::ihe4] += -0.5*yhe4*yhe4*yhe4*r3a;
405 dydt[Approx8Net::ihe4] += -yhe4*yo16*ro16a;
406 dydt[Approx8Net::ihe4] += 0.5*yc12*yc12*r1212;
407 dydt[Approx8Net::ihe4] += yc12*yo16*r1216;
408 dydt[Approx8Net::ihe4] += -yhe4*yne20*rne20a;
409
410 dydt[Approx8Net::ic12] = (1./6.)*yhe4*yhe4*yhe4*r3a;
411 dydt[Approx8Net::ic12] += -yhe4*yc12*rc12a;
412 dydt[Approx8Net::ic12] += -yh1*yc12*rc12p;
413 dydt[Approx8Net::ic12] += yh1*yn14*aFrac*rn14p;
414 dydt[Approx8Net::ic12] += -yc12*yc12*r1212;
415 dydt[Approx8Net::ic12] += -yc12*yo16*r1216;
416
417 dydt[Approx8Net::in14] = yh1*yc12*rc12p;
418 dydt[Approx8Net::in14] += -yh1*yn14*rn14p;
419 dydt[Approx8Net::in14] += yh1*yo16*ro16p;
420 dydt[Approx8Net::in14] += -yhe4*yn14*rn14a;
421
422 dydt[Approx8Net::io16] = yhe4*yc12*rc12a;
423 dydt[Approx8Net::io16] += yh1*yn14*pFrac*rn14p;
424 dydt[Approx8Net::io16] += -yh1*yo16*ro16p;
425 dydt[Approx8Net::io16] += -yc12*yo16*r1216;
426 dydt[Approx8Net::io16] += -yhe4*yo16*ro16a;
427
428 dydt[Approx8Net::ine20] = 0.5*yc12*yc12*r1212;
429 dydt[Approx8Net::ine20] += yhe4*yn14*rn14a;
430 dydt[Approx8Net::ine20] += yhe4*yo16*ro16a;
431 dydt[Approx8Net::ine20] += -yhe4*yne20*rne20a;
432
433 dydt[Approx8Net::img24] = yc12*yo16*r1216;
434 dydt[Approx8Net::img24] += yhe4*yne20*rne20a;
435
436 dydt[Approx8Net::iTemp] = 0.;
437 dydt[Approx8Net::iDensity] = 0.;
438
439 // energy accounting
440 double eNuc = 0.;
441 for (int i=0; i<Approx8Net::nIso; i++) {
442 eNuc += Approx8Net::mIon[i]*dydt[i];
443 }
444 dydt[Approx8Net::iEnergy] = -eNuc*avo*clight*clight;
445 }
446
448
450 m_y = convert_netIn(netIn);
451 m_tMax = netIn.tMax;
452 m_dt0 = netIn.dt0;
453
454 const double stiff_abs_tol = m_config.get<double>("Network:Approx8:Stiff:AbsTol", 1.0e-6);
455 const double stiff_rel_tol = m_config.get<double>("Network:Approx8:Stiff:RelTol", 1.0e-6);
456 const double nonstiff_abs_tol = m_config.get<double>("Network:Approx8:NonStiff:AbsTol", 1.0e-6);
457 const double nonstiff_rel_tol = m_config.get<double>("Network:Approx8:NonStiff:RelTol", 1.0e-6);
458
459 int num_steps = -1;
460
461 if (m_stiff) {
462 LOG_DEBUG(m_logger, "Using stiff solver for Approx8Network");
463 num_steps = integrate_const(
464 make_dense_output<rosenbrock4<double>>(stiff_abs_tol, stiff_rel_tol),
465 std::make_pair(ODE(), Jacobian()),
466 m_y,
467 0.0,
468 m_tMax,
469 m_dt0
470 );
471
472 } else {
473 LOG_DEBUG(m_logger, "Using non stiff solver for Approx8Network");
474 num_steps = integrate_const (
475 make_dense_output<runge_kutta_dopri5<vector_type>>(nonstiff_abs_tol, nonstiff_rel_tol),
476 ODE(),
477 m_y,
478 0.0,
479 m_tMax,
480 m_dt0
481 );
482 }
483
484 double ySum = 0.0;
485 for (int i = 0; i < Approx8Net::nIso; i++) {
486 m_y[i] *= Approx8Net::aIon[i];
487 ySum += m_y[i];
488 }
489 for (int i = 0; i < Approx8Net::nIso; i++) {
490 m_y[i] /= ySum;
491 }
492
493 NetOut netOut;
494 std::vector<double> outComposition;
495 outComposition.reserve(Approx8Net::nVar);
496
497 for (int i = 0; i < Approx8Net::nIso; i++) {
498 outComposition.push_back(m_y[i]);
499 }
500 netOut.energy = m_y[Approx8Net::iEnergy];
501 netOut.num_steps = num_steps;
502
503 const std::vector<std::string> symbols = {"H-1", "He-3", "He-4", "C-12", "N-14", "O-16", "Ne-20", "Mg-24"};
504 netOut.composition = serif::composition::Composition(symbols, outComposition);
505
506 return netOut;
507 }
508
509 void Approx8Network::setStiff(bool stiff) {
510 m_stiff = stiff;
511 }
512
515 y[Approx8Net::ih1] = netIn.composition.getMassFraction("H-1");
516 y[Approx8Net::ihe3] = netIn.composition.getMassFraction("He-3");
517 y[Approx8Net::ihe4] = netIn.composition.getMassFraction("He-4");
518 y[Approx8Net::ic12] = netIn.composition.getMassFraction("C-12");
519 y[Approx8Net::in14] = netIn.composition.getMassFraction("N-14");
520 y[Approx8Net::io16] = netIn.composition.getMassFraction("O-16");
521 y[Approx8Net::ine20] = netIn.composition.getMassFraction("Ne-20");
522 y[Approx8Net::img24] = netIn.composition.getMassFraction("Mg-24");
523 y[Approx8Net::iTemp] = netIn.temperature;
524 y[Approx8Net::iDensity] = netIn.density;
525 y[Approx8Net::iEnergy] = netIn.energy;
526
527 double ySum = 0.0;
528 for (int i = 0; i < Approx8Net::nIso; i++) {
529 y[i] /= Approx8Net::aIon[i];
530 ySum += y[i];
531 }
532 for (int i = 0; i < Approx8Net::nIso; i++) {
533 y[i] /= ySum;
534 }
535
536 return y;
537 }
538};
539
540
541// main program
542
Header file for the Approx8 nuclear reaction network.
Manages the composition of elements.
Class to manage a collection of constants.
Definition const.h:63
static Constants & getInstance()
get instance of constants singleton
Definition const.h:99
Constant get(const std::string &key) const
Get a constant by key.
Definition const.cpp:42
Network(const NetworkFormat format=NetworkFormat::APPROX8)
Definition network.cpp:32
quill::Logger * m_logger
Logger instance.
Definition network.h:136
serif::config::Config & m_config
Configuration instance.
Definition network.h:134
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
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
@ APPROX8
Approx8 nuclear reaction network format.
Definition network.h:37
const double value
Value of the constant.
Definition const.h:33
Input structure for the network evaluation.
Definition network.h:65
Output structure for the network evaluation.
Definition network.h:89
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