12 #include "sundials/sundials_types.h"
13 #include "sundials/sundials_math.h"
14 #include "sundials/sundials_nvector.h"
15 #include "nvector/nvector_serial.h"
16 #include "cvodes/cvodes.h"
17 #if CT_SUNDIALS_VERSION >= 30
18 #if CT_SUNDIALS_USE_LAPACK
19 #include "sunlinsol/sunlinsol_lapackdense.h"
20 #include "sunlinsol/sunlinsol_lapackband.h"
22 #include "sunlinsol/sunlinsol_dense.h"
23 #include "sunlinsol/sunlinsol_band.h"
25 #include "sunlinsol/sunlinsol_spgmr.h"
26 #include "cvodes/cvodes_direct.h"
27 #include "cvodes/cvodes_diag.h"
28 #include "cvodes/cvodes_spils.h"
30 #if CT_SUNDIALS_USE_LAPACK
31 #include "cvodes/cvodes_lapack.h"
33 #include "cvodes/cvodes_dense.h"
34 #include "cvodes/cvodes_band.h"
36 #include "cvodes/cvodes_diag.h"
37 #include "cvodes/cvodes_spgmr.h"
43 #if CT_SUNDIALS_VERSION < 25
44 typedef int sd_size_t;
46 typedef long int sd_size_t;
61 static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot,
void* f_data)
64 return f->
eval_nothrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
70 static void cvodes_err(
int error_code,
const char* module,
71 const char*
function,
char* msg,
void* eh_data)
79 CVodesIntegrator::CVodesIntegrator() :
101 m_maxErrTestFails(0),
104 m_mupper(0), m_mlower(0),
109 CVodesIntegrator::~CVodesIntegrator()
113 CVodeSensFree(m_cvode_mem);
115 CVodeFree(&m_cvode_mem);
118 #if CT_SUNDIALS_VERSION >= 30
119 SUNLinSolFree((SUNLinearSolver)
m_linsol);
124 N_VDestroy_Serial(m_y);
127 N_VDestroy_Serial(m_abstol);
130 N_VDestroy_Serial(m_dky);
133 N_VDestroyVectorArray_Serial(m_yS,
static_cast<sd_size_t
>(m_np));
139 return NV_Ith_S(m_y, k);
144 return NV_DATA_S(m_y);
153 N_VDestroy_Serial(m_abstol);
155 m_abstol = N_VNew_Serial(
static_cast<sd_size_t
>(n));
157 for (
size_t i=0; i<n; i++) {
158 NV_Ith_S(m_abstol, i) = abstol[i];
172 m_reltolsens = reltol;
173 m_abstolsens = abstol;
188 throw CanteraError(
"CVodesIntegrator::setMethod",
"unknown method");
196 CVodeSetMaxStep(m_cvode_mem, hmax);
204 CVodeSetMinStep(m_cvode_mem, hmin);
212 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
223 m_maxErrTestFails = n;
225 CVodeSetMaxErrTestFails(m_cvode_mem, n);
229 void CVodesIntegrator::sensInit(
double t0,
FuncEval& func)
234 N_Vector y = N_VNew_Serial(
static_cast<sd_size_t
>(func.
neq()));
235 m_yS = N_VCloneVectorArray_Serial(
static_cast<sd_size_t
>(m_np), y);
236 for (
size_t n = 0; n < m_np; n++) {
237 N_VConst(0.0, m_yS[n]);
239 N_VDestroy_Serial(y);
241 int flag = CVodeSensInit(m_cvode_mem,
static_cast<sd_size_t
>(m_np),
242 CV_STAGGERED, CVSensRhsFn(0), m_yS);
244 if (flag != CV_SUCCESS) {
245 throw CanteraError(
"CVodesIntegrator::sensInit",
"Error in CVodeSensInit");
248 for (
size_t n = 0; n < m_np; n++) {
253 flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
265 N_VDestroy_Serial(m_y);
267 m_y = N_VNew_Serial(
static_cast<sd_size_t
>(m_neq));
270 N_VDestroy_Serial(m_dky);
272 m_dky = N_VNew_Serial(
static_cast<sd_size_t
>(m_neq));
273 N_VConst(0.0, m_dky);
275 if (m_itol == CV_SV && m_nabs < m_neq) {
277 "not enough absolute tolerance values specified.");
283 CVodeFree(&m_cvode_mem);
289 #if CT_SUNDIALS_VERSION < 40
290 m_cvode_mem = CVodeCreate(m_method, CV_NEWTON);
292 m_cvode_mem = CVodeCreate(m_method);
296 "CVodeCreate failed.");
299 int flag = CVodeInit(m_cvode_mem,
cvodes_rhs, m_t0, m_y);
300 if (flag != CV_SUCCESS) {
301 if (flag == CV_MEM_FAIL) {
303 "Memory allocation failed.");
304 }
else if (flag == CV_ILL_INPUT) {
306 "Illegal value for CVodeInit input argument.");
309 "CVodeInit failed.");
312 CVodeSetErrHandlerFn(m_cvode_mem, &
cvodes_err,
this);
314 if (m_itol == CV_SV) {
315 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
317 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
319 if (flag != CV_SUCCESS) {
320 if (flag == CV_MEM_FAIL) {
322 "Memory allocation failed.");
323 }
else if (flag == CV_ILL_INPUT) {
325 "Illegal value for CVodeInit input argument.");
328 "CVodeInit failed.");
332 flag = CVodeSetUserData(m_cvode_mem, &func);
333 if (flag != CV_SUCCESS) {
335 "CVodeSetUserData failed.");
339 flag = CVodeSetSensParams(m_cvode_mem, func.
m_sens_params.data(),
341 if (flag != CV_SUCCESS) {
343 "CVodeSetSensParams failed.");
349 void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
357 int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
358 if (result != CV_SUCCESS) {
360 "CVodeReInit failed. result = {}", result);
367 if (m_type == DENSE + NOJAC) {
368 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
369 #if CT_SUNDIALS_VERSION >= 30
370 SUNLinSolFree((SUNLinearSolver)
m_linsol);
375 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
377 #if CT_SUNDIALS_USE_LAPACK
382 CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver)
m_linsol,
385 #if CT_SUNDIALS_USE_LAPACK
386 CVLapackDense(m_cvode_mem, N);
388 CVDense(m_cvode_mem, N);
391 }
else if (m_type == DIAG) {
393 }
else if (m_type == GMRES) {
394 #if CT_SUNDIALS_VERSION >= 30
395 m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
396 CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver)
m_linsol);
398 CVSpgmr(m_cvode_mem, PREC_NONE, 0);
400 }
else if (m_type == BAND + NOJAC) {
401 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
402 long int nu = m_mupper;
403 long int nl = m_mlower;
404 #if CT_SUNDIALS_VERSION >= 30
405 SUNLinSolFree((SUNLinearSolver)
m_linsol);
407 #if CT_SUNDIALS_VERSION < 40
414 "Unable to create SUNBandMatrix of size {} with bandwidths "
415 "{} and {}", N, nu, nl);
417 #if CT_SUNDIALS_USE_LAPACK
422 CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver)
m_linsol,
425 #if CT_SUNDIALS_USE_LAPACK
426 CVLapackBand(m_cvode_mem, N, nu, nl);
428 CVBand(m_cvode_mem, N, nu, nl);
433 "unsupported option");
437 CVodeSetMaxOrd(m_cvode_mem, m_maxord);
439 if (m_maxsteps > 0) {
440 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
443 CVodeSetMaxStep(m_cvode_mem, m_hmax);
446 CVodeSetMinStep(m_cvode_mem, m_hmin);
448 if (m_maxErrTestFails > 0) {
449 CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
458 int flag = CVode(m_cvode_mem, tout, m_y, &
m_time, CV_NORMAL);
459 if (flag != CV_SUCCESS) {
461 if (!f_errs.empty()) {
462 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
465 "CVodes error encountered. Error code: {}\n{}\n"
467 "Components with largest weighted error estimates:\n{}",
475 int flag = CVode(m_cvode_mem, tout, m_y, &
m_time, CV_ONE_STEP);
476 if (flag != CV_SUCCESS) {
478 if (!f_errs.empty()) {
479 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
482 "CVodes error encountered. Error code: {}\n{}\n"
484 "Components with largest weighted error estimates:\n{}",
494 int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
495 if (flag != CV_SUCCESS) {
497 if (!f_errs.empty()) {
498 f_errs =
"Exceptions caught evaluating derivative:\n" + f_errs;
501 "CVodes error encountered. Error code: {}\n{}\n"
505 return NV_DATA_S(m_dky);
511 CVodeGetLastOrder(m_cvode_mem, &ord);
518 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
522 double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
529 int flag = CVodeGetSens(m_cvode_mem, &
m_time, m_yS);
530 if (flag != CV_SUCCESS) {
531 throw CanteraError(
"CVodesIntegrator::sensitivity",
532 "CVodeGetSens failed. Error code: {}", flag);
538 throw CanteraError(
"CVodesIntegrator::sensitivity",
539 "sensitivity: k out of range ({})", k);
542 throw CanteraError(
"CVodesIntegrator::sensitivity",
543 "sensitivity: p out of range ({})", p);
545 return NV_Ith_S(m_yS[p],k);
550 N_Vector errs = N_VNew_Serial(
static_cast<sd_size_t
>(m_neq));
551 N_Vector errw = N_VNew_Serial(
static_cast<sd_size_t
>(m_neq));
552 CVodeGetErrWeights(m_cvode_mem, errw);
553 CVodeGetEstLocalErrors(m_cvode_mem, errs);
555 vector<tuple<double, double, size_t> > weightedErrors;
556 for (
size_t i=0; i<m_neq; i++) {
557 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
558 weightedErrors.emplace_back(-abs(err), err, i);
563 N = std::min(N,
static_cast<int>(m_neq));
564 sort(weightedErrors.begin(), weightedErrors.end());
565 fmt::memory_buffer s;
566 for (
int i=0; i<N; i++) {
567 format_to(s,
"{}: {}\n",
568 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));