12#include "cantera/numerics/sundials_headers.h"
18 return N_VNew_Serial(
static_cast<sd_size_t
>(N), context.get());
35 static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot,
void* f_data)
38 return f->
evalNoThrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
41 #if SUNDIALS_VERSION_MAJOR >= 7
46 static void sundials_err(
int line,
const char *func,
const char *file,
47 const char *msg, SUNErrCode err_code,
48 void *err_user_data, SUNContext sunctx)
50 CVodesIntegrator* integrator = (CVodesIntegrator*) err_user_data;
51 integrator->m_error_message = fmt::format(
"{}: {}\n", func, msg);
58 static void cvodes_err(
int error_code,
const char* module,
59 const char* function,
char* msg,
void* eh_data)
67 static int cvodes_prec_setup(sunrealtype t, N_Vector y, N_Vector ydot,
68 sunbooleantype jok, sunbooleantype *jcurPtr,
69 sunrealtype gamma,
void *f_data)
71 FuncEval* f = (FuncEval*) f_data;
74 return f->preconditioner_setup_nothrow(t, NV_DATA_S(y), gamma);
76 f->updatePreconditioner(gamma);
82 static int cvodes_prec_solve(sunrealtype t, N_Vector y, N_Vector ydot, N_Vector r,
83 N_Vector z, sunrealtype gamma, sunrealtype delta,
86 FuncEval* f = (FuncEval*) f_data;
87 return f->preconditioner_solve_nothrow(NV_DATA_S(r),NV_DATA_S(z));
99 static int cvodes_root(sunrealtype t, N_Vector y, sunrealtype* gout,
void* user_data)
101 auto* f =
static_cast<FuncEval*
>(user_data);
112CVodesIntegrator::~CVodesIntegrator()
121 SUNLinSolFree((SUNLinearSolver)
m_linsol);
125 N_VDestroy_Serial(
m_y);
128 N_VDestroy_Serial(m_abstol);
131 N_VDestroy_Serial(m_dky);
134 N_VDestroyVectorArray(m_yS,
static_cast<int>(m_np));
141 return NV_Ith_S(
m_y, k);
146 return NV_DATA_S(
m_y);
155 N_VDestroy_Serial(m_abstol);
159 for (
size_t i=0; i<n; i++) {
160 NV_Ith_S(m_abstol, i) = abstol[i];
164 int flag = CVodeSVtolerances(
m_cvode_mem, m_reltol, m_abstol);
165 checkError(flag,
"setTolerances",
"CVodeSVtolerances");
175 int flag = CVodeSStolerances(
m_cvode_mem, m_reltol, m_abstols);
176 checkError(flag,
"setTolerances",
"CVodeSStolerances");
182 m_reltolsens = reltol;
183 m_abstolsens = abstol;
185 vector<double> atol(m_np);
186 for (
size_t n = 0; n < m_np; n++) {
191 int flag = CVodeSensSStolerances(
m_cvode_mem, m_reltolsens, atol.data());
192 checkError(flag,
"setSensitivityTolerances",
"CVodeSensSStolerances");
203 throw CanteraError(
"CVodesIntegrator::setMethod",
"unknown method");
238 m_maxErrTestFails = n;
250 m_nRootFunctions = nroots;
256 CVRootFn root_cb = nroots ? &
cvodes_root :
nullptr;
257 int flag = CVodeRootInit(
m_cvode_mem,
static_cast<int>(nroots), root_cb);
258 checkError(flag,
"setRootFunctionCount",
"CVodeRootInit");
261void CVodesIntegrator::sensInit(
double t0,
FuncEval& func)
267 m_yS = N_VCloneVectorArray(
static_cast<int>(m_np), y);
268 for (
size_t n = 0; n < m_np; n++) {
269 N_VConst(0.0, m_yS[n]);
271 N_VDestroy_Serial(y);
273 int flag = CVodeSensInit(
m_cvode_mem,
static_cast<int>(m_np),
274 CV_STAGGERED, CVSensRhsFn(0), m_yS);
275 checkError(flag,
"sensInit",
"CVodeSensInit");
288 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
292 N_VDestroy_Serial(
m_y);
297 N_VDestroy_Serial(m_dky);
300 N_VConst(0.0, m_dky);
302 if (m_itol == CV_SV && m_nabs < m_neq) {
304 "not enough absolute tolerance values specified.");
316 "CVodeCreate failed.");
320 if (flag != CV_SUCCESS) {
321 if (flag == CV_MEM_FAIL) {
323 "Memory allocation failed.");
324 }
else if (flag == CV_ILL_INPUT) {
326 "Illegal value for CVodeInit input argument.");
329 "CVodeInit failed.");
332 #if SUNDIALS_VERSION_MAJOR >= 7
333 flag = SUNContext_PushErrHandler(
m_sundials_ctx.get(), &sundials_err,
this);
338 if (m_itol == CV_SV) {
339 flag = CVodeSVtolerances(
m_cvode_mem, m_reltol, m_abstol);
340 checkError(flag,
"initialize",
"CVodeSVtolerances");
342 flag = CVodeSStolerances(
m_cvode_mem, m_reltol, m_abstols);
343 checkError(flag,
"initialize",
"CVodeSStolerances");
347 checkError(flag,
"initialize",
"CVodeSetUserData");
349 m_nRootFunctions =
npos;
356 checkError(flag,
"initialize",
"CVodeSetSensParams");
361void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
370 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
374 checkError(result,
"reinitialize",
"CVodeReInit");
375 m_nRootFunctions =
npos;
382 if (m_type ==
"DENSE") {
383 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
384 SUNLinSolFree((SUNLinearSolver)
m_linsol);
389 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
392 #if CT_SUNDIALS_USE_LAPACK
403 "Error creating Sundials dense linear solver object");
404 }
else if (flag != CV_SUCCESS) {
406 "Error connecting linear solver to CVODES. "
407 "Sundials error code: {}", flag);
411 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
413 "Preconditioning is not available with the specified problem type.");
415 }
else if (m_type ==
"DIAG") {
418 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
420 "Preconditioning is not available with the specified problem type.");
422 }
else if (m_type ==
"GMRES") {
423 SUNLinSolFree((SUNLinearSolver)
m_linsol);
427 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
428 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver)
m_linsol,
430 CVodeSetPreconditioner(
m_cvode_mem, cvodes_prec_setup,
433 }
else if (m_type ==
"BAND") {
434 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
435 sd_size_t nu = m_mupper;
436 sd_size_t nl = m_mlower;
437 SUNLinSolFree((SUNLinearSolver)
m_linsol);
442 "Unable to create SUNBandMatrix of size {} with bandwidths "
443 "{} and {}", N, nu, nl);
445 #if CT_SUNDIALS_USE_LAPACK
456 "unsupported linear solver flag '{}'", m_type);
461 checkError(flag,
"applyOptions",
"CVodeSetMaxOrd");
463 if (m_maxsteps > 0) {
472 if (m_maxErrTestFails > 0) {
473 CVodeSetMaxErrTestFails(
m_cvode_mem, m_maxErrTestFails);
481 }
else if (tout <
m_time) {
483 "Cannot integrate backwards in time.\n"
484 "Requested time {} < current time {}",
489 if (nsteps >= m_maxsteps) {
491 if (!f_errs.empty()) {
492 f_errs =
"\nExceptions caught during RHS evaluation:\n" + f_errs;
495 "Maximum number of timesteps ({}) taken without reaching output "
496 "time ({}).\nCurrent integrator time: {}{}",
500 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
502 if (!f_errs.empty()) {
503 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
506 "CVodes error encountered. Error code: {}\n{}\n"
508 "Components with largest weighted error estimates:\n{}",
511 if (flag == CV_ROOT_RETURN) {
522 double t_eval = tout;
532 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
534 if (!f_errs.empty()) {
535 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
538 "CVodes error encountered. Error code: {}\n{}\n"
540 "Components with largest weighted error estimates:\n{}",
551 int flag = CVodeGetDky(
m_cvode_mem, tout, n, m_dky);
552 checkError(flag,
"derivative",
"CVodeGetDky");
553 return NV_DATA_S(m_dky);
576 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
577 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
578 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
579 nonlinConvFails = 0, orderReductions = 0, stepSolveFails = 0;
585 CVodeGetNonlinSolvStats(
m_cvode_mem, &nonlinIters, &nonlinConvFails);
586 CVodeGetNumErrTestFails(
m_cvode_mem, &errTestFails);
588 CVodeGetNumStabLimOrderReds(
m_cvode_mem, &orderReductions);
593 CVodeGetNumLinConvFails(
m_cvode_mem, &linConvFails);
596 CVodeGetNumJTSetupEvals(
m_cvode_mem, &jtSetupEvals);
598 CVodeGetNumStepSolveFails(
m_cvode_mem, &stepSolveFails);
600 stats[
"steps"] = steps;
601 stats[
"step_solve_fails"] = stepSolveFails;
602 stats[
"rhs_evals"] = rhsEvals;
603 stats[
"nonlinear_iters"] = nonlinIters;
604 stats[
"nonlinear_conv_fails"] = nonlinConvFails;
605 stats[
"err_test_fails"] = errTestFails;
607 stats[
"stab_order_reductions"] = orderReductions;
609 stats[
"jac_evals"] = jacEvals;
610 stats[
"lin_solve_setups"] = linSetup;
611 stats[
"lin_rhs_evals"] = linRhsEvals;
612 stats[
"lin_iters"] = linIters;
613 stats[
"lin_conv_fails"] = linConvFails;
614 stats[
"prec_evals"] = precEvals;
615 stats[
"prec_solves"] = precSolves;
616 stats[
"jt_vec_setup_evals"] = jtSetupEvals;
617 stats[
"jt_vec_prod_evals"] = jTimesEvals;
621double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
629 checkError(flag,
"sensitivity",
"CVodeGetSens");
634 throw CanteraError(
"CVodesIntegrator::sensitivity",
635 "sensitivity: k out of range ({})", k);
638 throw CanteraError(
"CVodesIntegrator::sensitivity",
639 "sensitivity: p out of range ({})", p);
641 return NV_Ith_S(m_yS[p],k);
651 vector<tuple<double, double, size_t>> weightedErrors;
652 for (
size_t i=0; i<m_neq; i++) {
653 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
654 weightedErrors.emplace_back(-abs(err), err, i);
659 N = std::min(N,
static_cast<int>(m_neq));
660 sort(weightedErrors.begin(), weightedErrors.end());
661 fmt::memory_buffer s;
662 for (
int i=0; i<N; i++) {
663 fmt_append(s,
"{}: {}\n",
664 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
670 const string& cvodesMethod)
const
672 if (flag == CV_SUCCESS) {
674 }
else if (flag == CV_MEM_NULL) {
676 "CVODES integrator is not initialized");
678 const char* flagname = CVodeGetReturnFlagName(flag);
680 "{} 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.
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.