GridFire 0.0.1a
General Purpose Nuclear Network
Loading...
Searching...
No Matches
engine_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 "fourdst/constants/const.h"
28#include "fourdst/config/config.h"
29#include "quill/LogMacros.h"
30
32#include "gridfire/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
75namespace gridfire::approx8{
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 double sum=0;
84 for (size_t i=0; i < a.size(); i++) {
85 sum += a[i] * b[i];
86 }
87 return sum;
88 }
89
90 // the fit to the nuclear reaction rates is of the form:
91 // exp( a0 + a1/T9 + a2/T9^(1/3) + a3*T9^(1/3) + a4*T9 + a5*T9^(5/3) + log(T9) )
92 // this function returns an array of the T9 terms in that order, where T9 is the temperatures in GigaKelvin
93 vec7 get_T9_array(const double &T) {
94 vec7 arr;
95 const double T9=1e-9*T;
96 const double T913=pow(T9,1./3.);
97
98 arr[0]=1;
99 arr[1]=1/T9;
100 arr[2]=1/T913;
101 arr[3]=T913;
102 arr[4]=T9;
103 arr[5]=pow(T9,5./3.);
104 arr[6]=log(T9);
105
106 return arr;
107 }
108
109 // this function uses the two preceding functions to evaluate the rate given the T9 array and the coefficients
110 double rate_fit(const vec7 &T9, const vec7 &coef){
111 return exp(sum_product(T9,coef));
112 }
113
114 // p + p -> d; this, like some of the other rates, this is a composite of multiple fits
115 double pp_rate(const vec7 &T9) {
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};
118 return rate_fit(T9,a1) + rate_fit(T9,a2);
119 }
120
121 // p + d -> he3
122 double dp_rate(const vec7 &T9) {
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};
125 return rate_fit(T9,a1) + rate_fit(T9,a2);
126 }
127
128 // he3 + he3 -> he4 + 2p
129 double he3he3_rate(const vec7 &T9){
130 constexpr vec7 a = {2.477880e+01,0,-12.27700,-0.1036990,-6.499670e-02,1.681910e-02,-6.666670e-01};
131 return rate_fit(T9,a);
132 }
133
134 // he3(he3,2p)he4 ** (missing both coefficients but have a reaction)
135 double he3he4_rate(const vec7 &T9){
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};
138 return rate_fit(T9,a1) + rate_fit(T9,a2);
139 }
140
141 // he4 + he4 + he4 -> c12 ** (missing middle coefficient but have other two)
142 double triple_alpha_rate(const vec7 &T9){
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};
146 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
147 }
148
149 // c12 + p -> n13
150 double c12p_rate(const vec7 &T9){
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};
153 return rate_fit(T9,a1) + rate_fit(T9,a2);
154 }
155
156 // c12 + he4 -> o16 ** (missing first coefficient but have the second)
157 double c12a_rate(const vec7 &T9){
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};
160 return rate_fit(T9,a1) + rate_fit(T9,a2);
161 }
162
163 // n14(p,g)o15 - o15 + p -> c12 + he4
164 double n14p_rate(const vec7 &T9){
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};
169 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
170 }
171
172 // n14(a,g)f18 assumed to go on to ne20
173 double n14a_rate(const vec7 &T9){
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};
177 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
178 }
179
180 // n15(p,a)c12 (CNO I)
181 double n15pa_rate(const vec7 &T9){
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};
186 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
187 }
188
189 // n15(p,g)o16 (CNO II)
190 double n15pg_rate(const vec7 &T9){
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};
194 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
195 }
196
197 double n15pg_frac(const vec7 &T9){
198 const double f1=n15pg_rate(T9);
199 const double f2=n15pa_rate(T9);
200 return f1/(f1+f2);
201 }
202
203 // o16(p,g)f17 then f17 -> o17(p,a)n14
204 double o16p_rate(const vec7 &T9){
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};
206 return rate_fit(T9,a);
207 }
208
209 // o16(a,g)ne20
210 double o16a_rate(const vec7 &T9){
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};
214 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3);
215 }
216
217 // ne20(a,g)mg24
218 double ne20a_rate(const vec7 &T9){
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};
223 return rate_fit(T9,a1) + rate_fit(T9,a2) + rate_fit(T9,a3) + rate_fit(T9,a4);
224 }
225
226 // c12(c12,a)ne20
227 double c12c12_rate(const vec7 &T9){
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};
229 return rate_fit(T9,a);
230 }
231
232 // c12(o16,a)mg24
233 double c12o16_rate(const vec7 &T9){
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};
235 return rate_fit(T9,a);
236 }
237
238
239 // for Boost ODE solvers either a struct or a class is required
240
241 // a Jacobian matrix for implicit solvers
242
243 void Jacobian::operator() ( const vector_type &y, matrix_type &J, double /* t */, vector_type &dfdt ) const {
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;
247 // EOS
249
250 // evaluate rates once per call
251 const double rpp=pp_rate(T9);
252 const double r33=he3he3_rate(T9);
253 const double r34=he3he4_rate(T9);
254 const double r3a=triple_alpha_rate(T9);
255 const double rc12p=c12p_rate(T9);
256 const double rc12a=c12a_rate(T9);
257 const double rn14p=n14p_rate(T9);
258 const double rn14a=n14a_rate(T9);
259 const double ro16p=o16p_rate(T9);
260 const double ro16a=o16a_rate(T9);
261 const double rne20a=ne20a_rate(T9);
262 const double r1212=c12c12_rate(T9);
263 const double r1216=c12o16_rate(T9);
264
265 const double pFrac=n15pg_frac(T9);
266 const double aFrac=1-pFrac;
267
268 const double yh1 = y[Approx8Net::ih1];
269 const double yhe3 = y[Approx8Net::ihe3];
270 const double yhe4 = y[Approx8Net::ihe4];
271 const double yc12 = y[Approx8Net::ic12];
272 const double yn14 = y[Approx8Net::in14];
273 const double yo16 = y[Approx8Net::io16];
274 const double yne20 = y[Approx8Net::ine20];
275
276 // zero all elements to begin
277 for (int i=0; i < Approx8Net::nVar; i++) {
278 for (int j=0; j < Approx8Net::nVar; j++) {
279 J(i,j)=0.0;
280 }
281 }
282
283 // h1 jacobian elements
284 J(Approx8Net::ih1,Approx8Net::ih1) = -3*yh1*rpp - 2*yc12*rc12p -2*yn14*rn14p -2*yo16*ro16p;
285 J(Approx8Net::ih1,Approx8Net::ihe3) = 2*yhe3*r33 - yhe4*r34;
286 J(Approx8Net::ih1,Approx8Net::ihe4) = -yhe3*r34;
287 J(Approx8Net::ih1,Approx8Net::ic12) = -2*yh1*rc12p;
288 J(Approx8Net::ih1,Approx8Net::in14) = -2*yh1*rn14p;
289 J(Approx8Net::ih1,Approx8Net::io16) = -2*yh1*ro16p;
290
291 // he3 jacobian elements
293 J(Approx8Net::ihe3,Approx8Net::ihe3) = -2*yhe3*r33 - yhe4*r34;
294 J(Approx8Net::ihe3,Approx8Net::ihe4) = -yhe3*r34;
295
296 // he4 jacobian elements
297 J(Approx8Net::ihe4,Approx8Net::ih1) = yn14*aFrac*rn14p + yo16*ro16p;
298 J(Approx8Net::ihe4,Approx8Net::ihe3) = yhe3*r33 - yhe4*r34;
299 J(Approx8Net::ihe4,Approx8Net::ihe4) = yhe3*r34 - 1.5*yhe4*yhe4*r3a - yc12*rc12a - 1.5*yn14*rn14a - yo16*ro16a - yne20*rne20a;
300 J(Approx8Net::ihe4,Approx8Net::ic12) = -yhe4*rc12a + yc12*r1212 + yo16*r1216;
301 J(Approx8Net::ihe4,Approx8Net::in14) = yh1*aFrac*rn14p - 1.5*yhe4*rn14a;
302 J(Approx8Net::ihe4,Approx8Net::io16) = yh1*ro16p - yhe4*ro16a + yc12*r1216;
303 J(Approx8Net::ihe4,Approx8Net::ine20) = -yhe4*rne20a;
304
305 // c12 jacobian elements
306 J(Approx8Net::ic12,Approx8Net::ih1) = -yc12*rc12p + yn14*aFrac*rn14p;
307 J(Approx8Net::ic12,Approx8Net::ihe4) = 0.5*yhe4*yhe4*r3a - yhe4*rc12a;
308 J(Approx8Net::ic12,Approx8Net::ic12) = -yh1*rc12p - yhe4*rc12a - yo16*r1216 - 2*yc12*r1212;
309 J(Approx8Net::ic12,Approx8Net::in14) = yh1*yn14*aFrac*rn14p;
310 J(Approx8Net::ic12,Approx8Net::io16) = -yc12*r1216;
311
312 // n14 jacobian elements
313 J(Approx8Net::in14,Approx8Net::ih1) = yc12*rc12p - yn14*rn14p + yo16*ro16p;
314 J(Approx8Net::in14,Approx8Net::ihe4) = -yn14*rn14a;
315 J(Approx8Net::in14,Approx8Net::ic12) = yh1*rc12p;
316 J(Approx8Net::in14,Approx8Net::in14) = -yh1*rn14p - yhe4*rn14a;
317 J(Approx8Net::in14,Approx8Net::io16) = yo16*ro16p;
318
319 // o16 jacobian elements
320 J(Approx8Net::io16,Approx8Net::ih1) = yn14*pFrac*rn14p - yo16*ro16p;
321 J(Approx8Net::io16,Approx8Net::ihe4) = yc12*rc12a - yo16*ro16a;
322 J(Approx8Net::io16,Approx8Net::ic12) = yhe4*rc12a - yo16*r1216;
323 J(Approx8Net::io16,Approx8Net::in14) = yh1*pFrac*rn14p;
324 J(Approx8Net::io16,Approx8Net::io16) = yh1*ro16p - yc12*r1216 -yhe4*ro16a;
325
326 // ne20 jacobian elements
327 J(Approx8Net::ine20,Approx8Net::ihe4) = yn14*rn14a + yo16*ro16a - yne20*rne20a;
328 J(Approx8Net::ine20,Approx8Net::ic12) = yc12*r1212;
329 J(Approx8Net::ine20,Approx8Net::in14) = yhe4*rn14a;
330 J(Approx8Net::ine20,Approx8Net::io16) = yo16*ro16a;
331 J(Approx8Net::ine20,Approx8Net::ine20) = -yhe4*rne20a;
332
333 // mg24 jacobian elements
334 J(Approx8Net::img24,Approx8Net::ihe4) = yne20*rne20a;
335 J(Approx8Net::img24,Approx8Net::ic12) = yo16*r1216;
336 J(Approx8Net::img24,Approx8Net::io16) = yc12*r1216;
337 J(Approx8Net::img24,Approx8Net::ine20) = yhe4*rne20a;
338
339 // energy accounting
340 for (int j=0; j<Approx8Net::nIso; j++) {
341 for (int i=0; i<Approx8Net::nIso; i++) {
342 J(Approx8Net::iEnergy,j) += J(i,j)*Approx8Net::mIon[i];
343 }
344 J(Approx8Net::iEnergy,j) *= -avo*clight*clight;
345 }
346 }
347
348 void ODE::operator() ( const vector_type &y, vector_type &dydt, double /* t */) const {
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;
352
353 // EOS
354 const double T = y[Approx8Net::iTemp];
355 const double den = y[Approx8Net::iDensity];
356 const vec7 T9=get_T9_array(T);
357
358 // rates
359 const double rpp=den*pp_rate(T9);
360 const double r33=den*he3he3_rate(T9);
361 const double r34=den*he3he4_rate(T9);
362 const double r3a=den*den*triple_alpha_rate(T9);
363 const double rc12p=den*c12p_rate(T9);
364 const double rc12a=den*c12a_rate(T9);
365 const double rn14p=den*n14p_rate(T9);
366 const double rn14a=n14a_rate(T9);
367 const double ro16p=den*o16p_rate(T9);
368 const double ro16a=den*o16a_rate(T9);
369 const double rne20a=den*ne20a_rate(T9);
370 const double r1212=den*c12c12_rate(T9);
371 const double r1216=den*c12o16_rate(T9);
372
373 const double pFrac=n15pg_frac(T9);
374 const double aFrac=1-pFrac;
375
376 const double yh1 = y[Approx8Net:: ih1];
377 const double yhe3 = y[Approx8Net:: ihe3];
378 const double yhe4 = y[Approx8Net:: ihe4];
379 const double yc12 = y[Approx8Net:: ic12];
380 const double yn14 = y[Approx8Net:: in14];
381 const double yo16 = y[Approx8Net:: io16];
382 const double yne20 = y[Approx8Net::ine20];
383
384 dydt[Approx8Net::ih1] = -1.5*yh1*yh1*rpp;
385 dydt[Approx8Net::ih1] += yhe3*yhe3*r33;
386 dydt[Approx8Net::ih1] += -yhe3*yhe4*r34;
387 dydt[Approx8Net::ih1] += -2*yh1*yc12*rc12p;
388 dydt[Approx8Net::ih1] += -2*yh1*yn14*rn14p;
389 dydt[Approx8Net::ih1] += -2*yh1*yo16*ro16p;
390
391 dydt[Approx8Net::ihe3] = 0.5*yh1*yh1*rpp;
392 dydt[Approx8Net::ihe3] += -yhe3*yhe3*r33;
393 dydt[Approx8Net::ihe3] += -yhe3*yhe4*r34;
394
395 dydt[Approx8Net::ihe4] = 0.5*yhe3*yhe3*r33;
396 dydt[Approx8Net::ihe4] += yhe3*yhe4*r34;
397 dydt[Approx8Net::ihe4] += -yhe4*yc12*rc12a;
398 dydt[Approx8Net::ihe4] += yh1*yn14*aFrac*rn14p;
399 dydt[Approx8Net::ihe4] += yh1*yo16*ro16p;
400 dydt[Approx8Net::ihe4] += -0.5*yhe4*yhe4*yhe4*r3a;
401 dydt[Approx8Net::ihe4] += -yhe4*yo16*ro16a;
402 dydt[Approx8Net::ihe4] += 0.5*yc12*yc12*r1212;
403 dydt[Approx8Net::ihe4] += yc12*yo16*r1216;
404 dydt[Approx8Net::ihe4] += -yhe4*yne20*rne20a;
405
406 dydt[Approx8Net::ic12] = (1./6.)*yhe4*yhe4*yhe4*r3a;
407 dydt[Approx8Net::ic12] += -yhe4*yc12*rc12a;
408 dydt[Approx8Net::ic12] += -yh1*yc12*rc12p;
409 dydt[Approx8Net::ic12] += yh1*yn14*aFrac*rn14p;
410 dydt[Approx8Net::ic12] += -yc12*yc12*r1212;
411 dydt[Approx8Net::ic12] += -yc12*yo16*r1216;
412
413 dydt[Approx8Net::in14] = yh1*yc12*rc12p;
414 dydt[Approx8Net::in14] += -yh1*yn14*rn14p;
415 dydt[Approx8Net::in14] += yh1*yo16*ro16p;
416 dydt[Approx8Net::in14] += -yhe4*yn14*rn14a;
417
418 dydt[Approx8Net::io16] = yhe4*yc12*rc12a;
419 dydt[Approx8Net::io16] += yh1*yn14*pFrac*rn14p;
420 dydt[Approx8Net::io16] += -yh1*yo16*ro16p;
421 dydt[Approx8Net::io16] += -yc12*yo16*r1216;
422 dydt[Approx8Net::io16] += -yhe4*yo16*ro16a;
423
424 dydt[Approx8Net::ine20] = 0.5*yc12*yc12*r1212;
425 dydt[Approx8Net::ine20] += yhe4*yn14*rn14a;
426 dydt[Approx8Net::ine20] += yhe4*yo16*ro16a;
427 dydt[Approx8Net::ine20] += -yhe4*yne20*rne20a;
428
429 dydt[Approx8Net::img24] = yc12*yo16*r1216;
430 dydt[Approx8Net::img24] += yhe4*yne20*rne20a;
431
432 dydt[Approx8Net::iTemp] = 0.;
433 dydt[Approx8Net::iDensity] = 0.;
434
435 // energy accounting
436 double eNuc = 0.;
437 for (int i=0; i<Approx8Net::nIso; i++) {
438 eNuc += Approx8Net::mIon[i]*dydt[i];
439 }
440 dydt[Approx8Net::iEnergy] = -eNuc*avo*clight*clight;
441 }
442
444
446 m_y = convert_netIn(netIn);
447 m_tMax = netIn.tMax;
448 m_dt0 = netIn.dt0;
449
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);
454
455 int num_steps = -1;
456
457 if (m_stiff) {
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),
461 std::make_pair(ODE(), Jacobian()),
462 m_y,
463 0.0,
464 m_tMax,
465 m_dt0
466 );
467
468 } else {
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),
472 ODE(),
473 m_y,
474 0.0,
475 m_tMax,
476 m_dt0
477 );
478 }
479
480 double ySum = 0.0;
481 for (int i = 0; i < Approx8Net::nIso; i++) {
482 m_y[i] *= Approx8Net::aIon[i];
483 ySum += m_y[i];
484 }
485 for (int i = 0; i < Approx8Net::nIso; i++) {
486 m_y[i] /= ySum;
487 }
488
489 NetOut netOut;
490 std::vector<double> outComposition;
491 outComposition.reserve(Approx8Net::nVar);
492
493 for (int i = 0; i < Approx8Net::nIso; i++) {
494 outComposition.push_back(m_y[i]);
495 }
497 netOut.num_steps = num_steps;
498
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);
501
502 return netOut;
503 }
504
505 void Approx8Network::setStiff(bool stiff) {
506 m_stiff = stiff;
507 }
508
511 y[Approx8Net::ih1] = netIn.composition.getNumberFraction("H-1");
512 y[Approx8Net::ihe3] = netIn.composition.getNumberFraction("He-3");
513 y[Approx8Net::ihe4] = netIn.composition.getNumberFraction("He-4");
514 y[Approx8Net::ic12] = netIn.composition.getNumberFraction("C-12");
515 y[Approx8Net::in14] = netIn.composition.getNumberFraction("N-14");
516 y[Approx8Net::io16] = netIn.composition.getNumberFraction("O-16");
517 y[Approx8Net::ine20] = netIn.composition.getNumberFraction("Ne-20");
518 y[Approx8Net::img24] = netIn.composition.getNumberFraction("Mg-24");
520 y[Approx8Net::iDensity] = netIn.density;
521 y[Approx8Net::iEnergy] = netIn.energy;
522
523 return y;
524 }
525};
526
527
528// main program
529
Network(const NetworkFormat format=NetworkFormat::APPROX8)
Definition network.cpp:41
quill::Logger * m_logger
Logger instance.
Definition network.h:98
fourdst::config::Config & m_config
Configuration instance.
Definition network.h:96
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.
Definition network.h:42
double density
Density in g/cm^3.
Definition network.h:58
double tMax
Maximum time.
Definition network.h:55
fourdst::composition::Composition composition
Composition of the network.
Definition network.h:54
double dt0
Initial time step.
Definition network.h:56
double temperature
Temperature in Kelvin.
Definition network.h:57
double energy
Energy in ergs.
Definition network.h:59
fourdst::composition::Composition composition
Composition of the network after evaluation.
Definition network.h:66
double energy
Energy in ergs after evaluation.
Definition network.h:68
int num_steps
Number of steps taken in the evaluation.
Definition network.h:67
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.