12#include "cantera/numerics/sundials_headers.h"
18#if SUNDIALS_VERSION_MAJOR >= 6
19 return N_VNew_Serial(
static_cast<sd_size_t
>(N), context.get());
21 return N_VNew_Serial(
static_cast<sd_size_t
>(N));
39 static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot,
void* f_data)
42 return f->
evalNoThrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
45 #if SUNDIALS_VERSION_MAJOR >= 7
50 static void sundials_err(
int line,
const char *func,
const char *file,
51 const char *msg, SUNErrCode err_code,
52 void *err_user_data, SUNContext sunctx)
54 CVodesIntegrator* integrator = (CVodesIntegrator*) err_user_data;
55 integrator->m_error_message = fmt::format(
"{}: {}\n", func, msg);
62 static void cvodes_err(
int error_code,
const char* module,
63 const char* function,
char* msg,
void* eh_data)
71 static int cvodes_prec_setup(sunrealtype t, N_Vector y, N_Vector ydot,
72 sunbooleantype jok, sunbooleantype *jcurPtr,
73 sunrealtype gamma,
void *f_data)
75 FuncEval* f = (FuncEval*) f_data;
78 return f->preconditioner_setup_nothrow(t, NV_DATA_S(y), gamma);
80 f->updatePreconditioner(gamma);
86 static int cvodes_prec_solve(sunrealtype t, N_Vector y, N_Vector ydot, N_Vector r,
87 N_Vector z, sunrealtype gamma, sunrealtype delta,
90 FuncEval* f = (FuncEval*) f_data;
91 return f->preconditioner_solve_nothrow(NV_DATA_S(r),NV_DATA_S(z));
103 static int cvodes_root(sunrealtype t, N_Vector y, sunrealtype* gout,
void* user_data)
105 auto* f =
static_cast<FuncEval*
>(user_data);
116CVodesIntegrator::~CVodesIntegrator()
125 SUNLinSolFree((SUNLinearSolver)
m_linsol);
129 N_VDestroy_Serial(
m_y);
132 N_VDestroy_Serial(m_abstol);
135 N_VDestroy_Serial(m_dky);
138 #if SUNDIALS_VERSION_MAJOR >= 6
139 N_VDestroyVectorArray(m_yS,
static_cast<int>(m_np));
141 N_VDestroyVectorArray_Serial(m_yS,
static_cast<int>(m_np));
149 return NV_Ith_S(
m_y, k);
154 return NV_DATA_S(
m_y);
163 N_VDestroy_Serial(m_abstol);
167 for (
size_t i=0; i<n; i++) {
168 NV_Ith_S(m_abstol, i) = abstol[i];
182 m_reltolsens = reltol;
183 m_abstolsens = abstol;
193 throw CanteraError(
"CVodesIntegrator::setMethod",
"unknown method");
228 m_maxErrTestFails = n;
240 m_nRootFunctions = nroots;
246 CVRootFn root_cb = nroots ? &
cvodes_root :
nullptr;
247 int flag = CVodeRootInit(
m_cvode_mem,
static_cast<int>(nroots), root_cb);
248 checkError(flag,
"setRootFunctionCount",
"CVodeRootInit");
251void CVodesIntegrator::sensInit(
double t0,
FuncEval& func)
257 #if SUNDIALS_VERSION_MAJOR >= 6
258 m_yS = N_VCloneVectorArray(
static_cast<int>(m_np), y);
260 m_yS = N_VCloneVectorArray_Serial(
static_cast<int>(m_np), y);
262 for (
size_t n = 0; n < m_np; n++) {
263 N_VConst(0.0, m_yS[n]);
265 N_VDestroy_Serial(y);
267 int flag = CVodeSensInit(
m_cvode_mem,
static_cast<int>(m_np),
268 CV_STAGGERED, CVSensRhsFn(0), m_yS);
269 checkError(flag,
"sensInit",
"CVodeSensInit");
271 vector<double> atol(m_np);
272 for (
size_t n = 0; n < m_np; n++) {
277 flag = CVodeSensSStolerances(
m_cvode_mem, m_reltolsens, atol.data());
278 checkError(flag,
"sensInit",
"CVodeSensSStolerances");
290 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
294 N_VDestroy_Serial(
m_y);
299 N_VDestroy_Serial(m_dky);
302 N_VConst(0.0, m_dky);
304 if (m_itol == CV_SV && m_nabs < m_neq) {
306 "not enough absolute tolerance values specified.");
318 #if SUNDIALS_VERSION_MAJOR < 6
325 "CVodeCreate failed.");
329 if (flag != CV_SUCCESS) {
330 if (flag == CV_MEM_FAIL) {
332 "Memory allocation failed.");
333 }
else if (flag == CV_ILL_INPUT) {
335 "Illegal value for CVodeInit input argument.");
338 "CVodeInit failed.");
341 #if SUNDIALS_VERSION_MAJOR >= 7
342 flag = SUNContext_PushErrHandler(
m_sundials_ctx.get(), &sundials_err,
this);
347 if (m_itol == CV_SV) {
348 flag = CVodeSVtolerances(
m_cvode_mem, m_reltol, m_abstol);
349 checkError(flag,
"initialize",
"CVodeSVtolerances");
351 flag = CVodeSStolerances(
m_cvode_mem, m_reltol, m_abstols);
352 checkError(flag,
"initialize",
"CVodeSStolerances");
356 checkError(flag,
"initialize",
"CVodeSetUserData");
358 m_nRootFunctions =
npos;
365 checkError(flag,
"initialize",
"CVodeSetSensParams");
370void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
379 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
383 checkError(result,
"reinitialize",
"CVodeReInit");
384 m_nRootFunctions =
npos;
391 if (m_type ==
"DENSE") {
392 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
393 SUNLinSolFree((SUNLinearSolver)
m_linsol);
395 #if SUNDIALS_VERSION_MAJOR >= 6
402 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
405 #if SUNDIALS_VERSION_MAJOR >= 6
406 #if CT_SUNDIALS_USE_LAPACK
416 #if CT_SUNDIALS_USE_LAPACK
426 "Error creating Sundials dense linear solver object");
427 }
else if (flag != CV_SUCCESS) {
429 "Error connecting linear solver to CVODES. "
430 "Sundials error code: {}", flag);
434 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
436 "Preconditioning is not available with the specified problem type.");
438 }
else if (m_type ==
"DIAG") {
441 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
443 "Preconditioning is not available with the specified problem type.");
445 }
else if (m_type ==
"GMRES") {
446 SUNLinSolFree((SUNLinearSolver)
m_linsol);
447 #if SUNDIALS_VERSION_MAJOR >= 6
455 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
456 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver)
m_linsol,
458 CVodeSetPreconditioner(
m_cvode_mem, cvodes_prec_setup,
461 }
else if (m_type ==
"BAND") {
462 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
463 sd_size_t nu = m_mupper;
464 sd_size_t nl = m_mlower;
465 SUNLinSolFree((SUNLinearSolver)
m_linsol);
467 #if SUNDIALS_VERSION_MAJOR >= 6
474 "Unable to create SUNBandMatrix of size {} with bandwidths "
475 "{} and {}", N, nu, nl);
477 #if SUNDIALS_VERSION_MAJOR >= 6
478 #if CT_SUNDIALS_USE_LAPACK
488 #if CT_SUNDIALS_USE_LAPACK
498 "unsupported linear solver flag '{}'", m_type);
503 checkError(flag,
"applyOptions",
"CVodeSetMaxOrd");
505 if (m_maxsteps > 0) {
514 if (m_maxErrTestFails > 0) {
515 CVodeSetMaxErrTestFails(
m_cvode_mem, m_maxErrTestFails);
523 }
else if (tout <
m_time) {
525 "Cannot integrate backwards in time.\n"
526 "Requested time {} < current time {}",
531 if (nsteps >= m_maxsteps) {
533 if (!f_errs.empty()) {
534 f_errs =
"\nExceptions caught during RHS evaluation:\n" + f_errs;
537 "Maximum number of timesteps ({}) taken without reaching output "
538 "time ({}).\nCurrent integrator time: {}{}",
542 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
544 if (!f_errs.empty()) {
545 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
548 "CVodes error encountered. Error code: {}\n{}\n"
550 "Components with largest weighted error estimates:\n{}",
553 if (flag == CV_ROOT_RETURN) {
564 double t_eval = tout;
574 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
576 if (!f_errs.empty()) {
577 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
580 "CVodes error encountered. Error code: {}\n{}\n"
582 "Components with largest weighted error estimates:\n{}",
593 int flag = CVodeGetDky(
m_cvode_mem, tout, n, m_dky);
594 checkError(flag,
"derivative",
"CVodeGetDky");
595 return NV_DATA_S(m_dky);
618 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
619 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
620 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
621 nonlinConvFails = 0, orderReductions = 0;
627 CVodeGetNonlinSolvStats(
m_cvode_mem, &nonlinIters, &nonlinConvFails);
628 CVodeGetNumErrTestFails(
m_cvode_mem, &errTestFails);
630 CVodeGetNumStabLimOrderReds(
m_cvode_mem, &orderReductions);
635 CVodeGetNumLinConvFails(
m_cvode_mem, &linConvFails);
638 CVodeGetNumJTSetupEvals(
m_cvode_mem, &jtSetupEvals);
641 #if SUNDIALS_VERSION_MAJOR >= 7 || (SUNDIALS_VERSION_MAJOR == 6 && SUNDIALS_VERSION_MINOR >= 2)
642 long int stepSolveFails = 0;
643 CVodeGetNumStepSolveFails(
m_cvode_mem, &stepSolveFails);
644 stats[
"step_solve_fails"] = stepSolveFails;
647 stats[
"steps"] = steps;
648 stats[
"rhs_evals"] = rhsEvals;
649 stats[
"nonlinear_iters"] = nonlinIters;
650 stats[
"nonlinear_conv_fails"] = nonlinConvFails;
651 stats[
"err_test_fails"] = errTestFails;
653 stats[
"stab_order_reductions"] = orderReductions;
655 stats[
"jac_evals"] = jacEvals;
656 stats[
"lin_solve_setups"] = linSetup;
657 stats[
"lin_rhs_evals"] = linRhsEvals;
658 stats[
"lin_iters"] = linIters;
659 stats[
"lin_conv_fails"] = linConvFails;
660 stats[
"prec_evals"] = precEvals;
661 stats[
"prec_solves"] = precSolves;
662 stats[
"jt_vec_setup_evals"] = jtSetupEvals;
663 stats[
"jt_vec_prod_evals"] = jTimesEvals;
667double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
675 checkError(flag,
"sensitivity",
"CVodeGetSens");
680 throw CanteraError(
"CVodesIntegrator::sensitivity",
681 "sensitivity: k out of range ({})", k);
684 throw CanteraError(
"CVodesIntegrator::sensitivity",
685 "sensitivity: p out of range ({})", p);
687 return NV_Ith_S(m_yS[p],k);
697 vector<tuple<double, double, size_t>> weightedErrors;
698 for (
size_t i=0; i<m_neq; i++) {
699 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
700 weightedErrors.emplace_back(-abs(err), err, i);
705 N = std::min(N,
static_cast<int>(m_neq));
706 sort(weightedErrors.begin(), weightedErrors.end());
707 fmt::memory_buffer s;
708 for (
int i=0; i<N; i++) {
709 fmt_append(s,
"{}: {}\n",
710 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
716 const string& cvodesMethod)
const
718 if (flag == CV_SUCCESS) {
720 }
else if (flag == CV_MEM_NULL) {
722 "CVODES integrator is not initialized");
724 const char* flagname = CVodeGetReturnFlagName(flag);
726 "{} returned error code {} ({}):\n{}",
A map of string keys to values whose type can vary at runtime.
Wrapper class for 'cvodes' integrator from LLNL.
double step(double tout) override
Integrate the system of equations.
void setMaxStepSize(double hmax) override
Set the maximum step size.
double * solution() override
The current value of the solution of the system of equations.
N_Vector m_y
The system state at m_time.
void * m_cvode_mem
CVODES internal memory object.
void setRootFunctionCount(size_t nroots) override
Configure how many event/root functions the integrator should monitor.
void * m_linsol
Sundials linear solver object.
int nEvals() const override
The number of function evaluations.
CVodesIntegrator()
Constructor.
double m_time
The current system time, corresponding to m_y.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS have been updated for at the current integrator ti...
int maxSteps() override
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
SundialsContext m_sundials_ctx
SUNContext object for Sundials>=6.0.
void checkError(long flag, const string &ctMethod, const string &cvodesMethod) const
Check whether a CVODES method indicated an error.
void setSensitivityTolerances(double reltol, double abstol) override
Set the sensitivity error tolerances.
void applyOptions()
Applies user-specified options to the underlying CVODES solver.
int lastOrder() const override
Order used during the last solution step.
void setMinStepSize(double hmin) override
Set the minimum step size.
void integrate(double tout) override
Integrate the system of equations.
void setTolerances(double reltol, size_t n, double *abstol) override
Set error tolerances.
double * derivative(double tout, int n) override
n-th derivative of the output function at time tout.
string m_error_message
Error message information provide by CVodes.
void setMaxErrTestFails(int n) override
Set the maximum permissible number of error test failures.
void setMaxSteps(int nmax) override
Set the maximum number of time-steps the integrator can take before reaching the next output time.
AnyMap solverStats() const override
Get solver stats from integrator.
void * m_linsol_matrix
matrix used by Sundials
double m_tInteg
The latest time reached by the integrator. May be greater than m_time.
void initialize(double t0, FuncEval &func) override
Initialize the integrator for a new problem.
string getErrorInfo(int N)
Returns a string listing the weighted error estimates associated with each solution component.
void setMethod(MethodType t) override
Set the solution method.
Base class for exceptions thrown by Cantera classes.
Virtual base class for ODE/DAE right-hand-side function evaluators.
int evalRootFunctionsNoThrow(double t, const double *y, double *gout)
Wrapper for evalRootFunctions that converts exceptions to return codes.
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
void clearErrors()
Clear any previously-stored suppressed errors.
int evalNoThrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
virtual size_t neq() const =0
Number of equations.
virtual size_t nRootFunctions() const
Number of event/root functions exposed to the integrator.
virtual size_t nparams() const
Number of sensitivity parameters.
string getErrors() const
Return a string containing the text of any suppressed errors.
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
shared_ptr< SystemJacobian > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
A wrapper for managing a SUNContext object, need for Sundials >= 6.0.
static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
MethodType
Specifies the method used to integrate the system of equations.
@ BDF_Method
Backward Differentiation.
static int cvodes_root(sunrealtype t, N_Vector y, sunrealtype *gout, void *user_data)
SUNDIALS callback that forwards root evaluations to the FuncEval.
static void cvodes_err(int error_code, const char *module, const char *function, char *msg, void *eh_data)
Function called by CVodes when an error is encountered instead of writing to stdout.
Contains declarations for string manipulation functions within Cantera.