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));
49 static void cvodes_err(
int error_code,
const char* module,
50 const char* function,
char* msg,
void* eh_data)
61 #if SUNDIALS_VERSION_MAJOR >= 7
62 static void sundials_err(
int line,
const char *func,
const char *file,
63 const char *msg, SUNErrCode err_code,
64 void *err_user_data, SUNContext sunctx)
66 CVodesIntegrator* integrator = (CVodesIntegrator*) err_user_data;
67 integrator->m_error_message = fmt::format(
"{}: {}\n", func, msg);
72 sunbooleantype jok, sunbooleantype *jcurPtr,
73 sunrealtype gamma,
void *f_data)
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));
101CVodesIntegrator::~CVodesIntegrator()
110 SUNLinSolFree((SUNLinearSolver)
m_linsol);
114 N_VDestroy_Serial(
m_y);
117 N_VDestroy_Serial(m_abstol);
120 N_VDestroy_Serial(m_dky);
123 #if SUNDIALS_VERSION_MAJOR >= 6
124 N_VDestroyVectorArray(m_yS,
static_cast<int>(m_np));
126 N_VDestroyVectorArray_Serial(m_yS,
static_cast<int>(m_np));
134 return NV_Ith_S(
m_y, k);
139 return NV_DATA_S(
m_y);
148 N_VDestroy_Serial(m_abstol);
152 for (
size_t i=0; i<n; i++) {
153 NV_Ith_S(m_abstol, i) = abstol[i];
167 m_reltolsens = reltol;
168 m_abstolsens = abstol;
178 throw CanteraError(
"CVodesIntegrator::setMethod",
"unknown method");
213 m_maxErrTestFails = n;
219void CVodesIntegrator::sensInit(
double t0,
FuncEval& func)
225 #if SUNDIALS_VERSION_MAJOR >= 6
226 m_yS = N_VCloneVectorArray(
static_cast<int>(m_np), y);
228 m_yS = N_VCloneVectorArray_Serial(
static_cast<int>(m_np), y);
230 for (
size_t n = 0; n < m_np; n++) {
231 N_VConst(0.0, m_yS[n]);
233 N_VDestroy_Serial(y);
235 int flag = CVodeSensInit(
m_cvode_mem,
static_cast<int>(m_np),
236 CV_STAGGERED, CVSensRhsFn(0), m_yS);
237 checkError(flag,
"sensInit",
"CVodeSensInit");
239 vector<double> atol(m_np);
240 for (
size_t n = 0; n < m_np; n++) {
245 flag = CVodeSensSStolerances(
m_cvode_mem, m_reltolsens, atol.data());
246 checkError(flag,
"sensInit",
"CVodeSensSStolerances");
258 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
262 N_VDestroy_Serial(
m_y);
267 N_VDestroy_Serial(m_dky);
270 N_VConst(0.0, m_dky);
272 if (m_itol == CV_SV && m_nabs < m_neq) {
274 "not enough absolute tolerance values specified.");
286 #if SUNDIALS_VERSION_MAJOR < 6
293 "CVodeCreate failed.");
297 if (flag != CV_SUCCESS) {
298 if (flag == CV_MEM_FAIL) {
300 "Memory allocation failed.");
301 }
else if (flag == CV_ILL_INPUT) {
303 "Illegal value for CVodeInit input argument.");
306 "CVodeInit failed.");
309 #if SUNDIALS_VERSION_MAJOR >= 7
310 flag = SUNContext_PushErrHandler(
m_sundials_ctx.get(), &sundials_err,
this);
315 if (m_itol == CV_SV) {
316 flag = CVodeSVtolerances(
m_cvode_mem, m_reltol, m_abstol);
317 checkError(flag,
"initialize",
"CVodeSVtolerances");
319 flag = CVodeSStolerances(
m_cvode_mem, m_reltol, m_abstols);
320 checkError(flag,
"initialize",
"CVodeSStolerances");
324 checkError(flag,
"initialize",
"CVodeSetUserData");
330 checkError(flag,
"initialize",
"CVodeSetSensParams");
335void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
344 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
348 checkError(result,
"reinitialize",
"CVodeReInit");
354 if (m_type ==
"DENSE") {
355 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
356 SUNLinSolFree((SUNLinearSolver)
m_linsol);
358 #if SUNDIALS_VERSION_MAJOR >= 6
365 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
368 #if SUNDIALS_VERSION_MAJOR >= 6
369 #if CT_SUNDIALS_USE_LAPACK
379 #if CT_SUNDIALS_USE_LAPACK
389 "Error creating Sundials dense linear solver object");
390 }
else if (flag != CV_SUCCESS) {
392 "Error connecting linear solver to CVODES. "
393 "Sundials error code: {}", flag);
397 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
399 "Preconditioning is not available with the specified problem type.");
401 }
else if (m_type ==
"DIAG") {
404 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
406 "Preconditioning is not available with the specified problem type.");
408 }
else if (m_type ==
"GMRES") {
409 SUNLinSolFree((SUNLinearSolver)
m_linsol);
410 #if SUNDIALS_VERSION_MAJOR >= 6
418 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
419 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver)
m_linsol,
424 }
else if (m_type ==
"BAND") {
425 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
426 sd_size_t nu = m_mupper;
427 sd_size_t nl = m_mlower;
428 SUNLinSolFree((SUNLinearSolver)
m_linsol);
430 #if SUNDIALS_VERSION_MAJOR >= 6
437 "Unable to create SUNBandMatrix of size {} with bandwidths "
438 "{} and {}", N, nu, nl);
440 #if SUNDIALS_VERSION_MAJOR >= 6
441 #if CT_SUNDIALS_USE_LAPACK
451 #if CT_SUNDIALS_USE_LAPACK
461 "unsupported linear solver flag '{}'", m_type);
466 checkError(flag,
"applyOptions",
"CVodeSetMaxOrd");
468 if (m_maxsteps > 0) {
477 if (m_maxErrTestFails > 0) {
478 CVodeSetMaxErrTestFails(
m_cvode_mem, m_maxErrTestFails);
486 }
else if (tout <
m_time) {
488 "Cannot integrate backwards in time.\n"
489 "Requested time {} < current time {}",
494 if (nsteps >= m_maxsteps) {
496 if (!f_errs.empty()) {
497 f_errs =
"\nExceptions caught during RHS evaluation:\n" + f_errs;
500 "Maximum number of timesteps ({}) taken without reaching output "
501 "time ({}).\nCurrent integrator time: {}{}",
505 if (flag != CV_SUCCESS) {
507 if (!f_errs.empty()) {
508 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
511 "CVodes error encountered. Error code: {}\n{}\n"
513 "Components with largest weighted error estimates:\n{}",
527 if (flag != CV_SUCCESS) {
529 if (!f_errs.empty()) {
530 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
533 "CVodes error encountered. Error code: {}\n{}\n"
535 "Components with largest weighted error estimates:\n{}",
546 int flag = CVodeGetDky(
m_cvode_mem, tout, n, m_dky);
547 checkError(flag,
"derivative",
"CVodeGetDky");
548 return NV_DATA_S(m_dky);
571 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
572 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
573 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
574 nonlinConvFails = 0, orderReductions = 0;
580 CVodeGetNonlinSolvStats(
m_cvode_mem, &nonlinIters, &nonlinConvFails);
581 CVodeGetNumErrTestFails(
m_cvode_mem, &errTestFails);
583 CVodeGetNumStabLimOrderReds(
m_cvode_mem, &orderReductions);
588 CVodeGetNumLinConvFails(
m_cvode_mem, &linConvFails);
591 CVodeGetNumJTSetupEvals(
m_cvode_mem, &jtSetupEvals);
594 #if SUNDIALS_VERSION_MAJOR >= 7 || (SUNDIALS_VERSION_MAJOR == 6 && SUNDIALS_VERSION_MINOR >= 2)
595 long int stepSolveFails = 0;
596 CVodeGetNumStepSolveFails(
m_cvode_mem, &stepSolveFails);
597 stats[
"step_solve_fails"] = stepSolveFails;
600 stats[
"steps"] = steps;
601 stats[
"rhs_evals"] = rhsEvals;
602 stats[
"nonlinear_iters"] = nonlinIters;
603 stats[
"nonlinear_conv_fails"] = nonlinConvFails;
604 stats[
"err_test_fails"] = errTestFails;
606 stats[
"stab_order_reductions"] = orderReductions;
608 stats[
"jac_evals"] = jacEvals;
609 stats[
"lin_solve_setups"] = linSetup;
610 stats[
"lin_rhs_evals"] = linRhsEvals;
611 stats[
"lin_iters"] = linIters;
612 stats[
"lin_conv_fails"] = linConvFails;
613 stats[
"prec_evals"] = precEvals;
614 stats[
"prec_solves"] = precSolves;
615 stats[
"jt_vec_setup_evals"] = jtSetupEvals;
616 stats[
"jt_vec_prod_evals"] = jTimesEvals;
620double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
628 checkError(flag,
"sensitivity",
"CVodeGetSens");
633 throw CanteraError(
"CVodesIntegrator::sensitivity",
634 "sensitivity: k out of range ({})", k);
637 throw CanteraError(
"CVodesIntegrator::sensitivity",
638 "sensitivity: p out of range ({})", p);
640 return NV_Ith_S(m_yS[p],k);
650 vector<tuple<double, double, size_t>> weightedErrors;
651 for (
size_t i=0; i<m_neq; i++) {
652 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
653 weightedErrors.emplace_back(-abs(err), err, i);
658 N = std::min(N,
static_cast<int>(m_neq));
659 sort(weightedErrors.begin(), weightedErrors.end());
660 fmt::memory_buffer s;
661 for (
int i=0; i<N; i++) {
662 fmt_append(s,
"{}: {}\n",
663 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
669 const string& cvodesMethod)
const
671 if (flag == CV_SUCCESS) {
673 }
else if (flag == CV_MEM_NULL) {
675 "CVODES integrator is not initialized");
677 const char* flagname = CVodeGetReturnFlagName(flag);
679 "{} 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 * 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.
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
int preconditioner_setup_nothrow(double t, double *y, double gamma)
Preconditioner setup that doesn't throw an error but returns a CVODES flag.
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 void updatePreconditioner(double gamma)
Update the preconditioner based on already computed jacobian values.
virtual size_t neq() const =0
Number of equations.
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.
static int cvodes_prec_setup(sunrealtype t, N_Vector y, N_Vector ydot, sunbooleantype jok, sunbooleantype *jcurPtr, sunrealtype gamma, void *f_data)
Function called by CVodes when an error is encountered instead of writing to stdout.
MethodType
Specifies the method used to integrate the system of equations.
@ BDF_Method
Backward Differentiation.
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.