12#include "cantera/numerics/sundials_headers.h"
18#if CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION >= 70
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 CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION < 40
288 #elif CT_SUNDIALS_VERSION < 60
295 "CVodeCreate failed.");
299 if (flag != CV_SUCCESS) {
300 if (flag == CV_MEM_FAIL) {
302 "Memory allocation failed.");
303 }
else if (flag == CV_ILL_INPUT) {
305 "Illegal value for CVodeInit input argument.");
308 "CVodeInit failed.");
311 #if CT_SUNDIALS_VERSION >= 70
312 flag = SUNContext_PushErrHandler(
m_sundials_ctx.get(), &sundials_err,
this);
317 if (m_itol == CV_SV) {
318 flag = CVodeSVtolerances(
m_cvode_mem, m_reltol, m_abstol);
319 checkError(flag,
"initialize",
"CVodeSVtolerances");
321 flag = CVodeSStolerances(
m_cvode_mem, m_reltol, m_abstols);
322 checkError(flag,
"initialize",
"CVodeSStolerances");
326 checkError(flag,
"initialize",
"CVodeSetUserData");
332 checkError(flag,
"initialize",
"CVodeSetSensParams");
337void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
346 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
350 checkError(result,
"reinitialize",
"CVodeReInit");
356 if (m_type ==
"DENSE") {
357 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
358 SUNLinSolFree((SUNLinearSolver)
m_linsol);
360 #if CT_SUNDIALS_VERSION >= 60
367 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
370 #if CT_SUNDIALS_VERSION >= 60
371 #if CT_SUNDIALS_USE_LAPACK
381 #if CT_SUNDIALS_USE_LAPACK
391 "Error creating Sundials dense linear solver object");
392 }
else if (flag != CV_SUCCESS) {
394 "Error connecting linear solver to CVODES. "
395 "Sundials error code: {}", flag);
399 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
401 "Preconditioning is not available with the specified problem type.");
403 }
else if (m_type ==
"DIAG") {
406 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
408 "Preconditioning is not available with the specified problem type.");
410 }
else if (m_type ==
"GMRES") {
411 #if CT_SUNDIALS_VERSION >= 60
414 #elif CT_SUNDIALS_VERSION >= 40
422 #if CT_SUNDIALS_VERSION >= 40
423 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
424 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver)
m_linsol,
430 if (
m_prec_side != PreconditionerSide::NO_PRECONDITION) {
431 SUNSPGMRSetPrecType((SUNLinearSolver)
m_linsol,
437 }
else if (m_type ==
"BAND") {
438 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
439 sd_size_t nu = m_mupper;
440 sd_size_t nl = m_mlower;
441 SUNLinSolFree((SUNLinearSolver)
m_linsol);
443 #if CT_SUNDIALS_VERSION >= 60
445 #elif CT_SUNDIALS_VERSION >= 40
452 "Unable to create SUNBandMatrix of size {} with bandwidths "
453 "{} and {}", N, nu, nl);
455 #if CT_SUNDIALS_VERSION >= 60
456 #if CT_SUNDIALS_USE_LAPACK
466 #if CT_SUNDIALS_USE_LAPACK
476 "unsupported linear solver flag '{}'", m_type);
481 checkError(flag,
"applyOptions",
"CVodeSetMaxOrd");
483 if (m_maxsteps > 0) {
492 if (m_maxErrTestFails > 0) {
493 CVodeSetMaxErrTestFails(
m_cvode_mem, m_maxErrTestFails);
501 }
else if (tout <
m_time) {
503 "Cannot integrate backwards in time.\n"
504 "Requested time {} < current time {}",
509 if (nsteps >= m_maxsteps) {
511 "Maximum number of timesteps ({}) taken without reaching output "
512 "time ({}).\nCurrent integrator time: {}",
516 if (flag != CV_SUCCESS) {
518 if (!f_errs.empty()) {
519 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
522 "CVodes error encountered. Error code: {}\n{}\n"
524 "Components with largest weighted error estimates:\n{}",
538 if (flag != CV_SUCCESS) {
540 if (!f_errs.empty()) {
541 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
544 "CVodes error encountered. Error code: {}\n{}\n"
546 "Components with largest weighted error estimates:\n{}",
557 int flag = CVodeGetDky(
m_cvode_mem, tout, n, m_dky);
558 checkError(flag,
"derivative",
"CVodeGetDky");
559 return NV_DATA_S(m_dky);
582 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
583 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
584 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
585 nonlinConvFails = 0, orderReductions = 0;
588 #if CT_SUNDIALS_VERSION >= 40
591 CVodeGetNonlinSolvStats(
m_cvode_mem, &nonlinIters, &nonlinConvFails);
592 CVodeGetNumErrTestFails(
m_cvode_mem, &errTestFails);
594 CVodeGetNumStabLimOrderReds(
m_cvode_mem, &orderReductions);
599 CVodeGetNumLinConvFails(
m_cvode_mem, &linConvFails);
602 CVodeGetNumJTSetupEvals(
m_cvode_mem, &jtSetupEvals);
605 warn_user(
"CVodesIntegrator::solverStats",
"Function not"
606 "supported with sundials versions less than 4.");
609 #if CT_SUNDIALS_VERION >= 60
610 long int stepSolveFails = 0;
611 CVodeGetNumStepSolveFails(
m_cvode_mem, &stepSolveFails);
612 stats[
"step_solve_fails"] = stepSolveFails;
615 stats[
"steps"] = steps;
616 stats[
"rhs_evals"] = rhsEvals;
617 stats[
"nonlinear_iters"] = nonlinIters;
618 stats[
"nonlinear_conv_fails"] = nonlinConvFails;
619 stats[
"err_test_fails"] = errTestFails;
621 stats[
"stab_order_reductions"] = orderReductions;
623 stats[
"jac_evals"] = jacEvals;
624 stats[
"lin_solve_setups"] = linSetup;
625 stats[
"lin_rhs_evals"] = linRhsEvals;
626 stats[
"lin_iters"] = linIters;
627 stats[
"lin_conv_fails"] = linConvFails;
628 stats[
"prec_evals"] = precEvals;
629 stats[
"prec_solves"] = precSolves;
630 stats[
"jt_vec_setup_evals"] = jtSetupEvals;
631 stats[
"jt_vec_prod_evals"] = jTimesEvals;
635double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
643 checkError(flag,
"sensitivity",
"CVodeGetSens");
648 throw CanteraError(
"CVodesIntegrator::sensitivity",
649 "sensitivity: k out of range ({})", k);
652 throw CanteraError(
"CVodesIntegrator::sensitivity",
653 "sensitivity: p out of range ({})", p);
655 return NV_Ith_S(m_yS[p],k);
665 vector<tuple<double, double, size_t>> weightedErrors;
666 for (
size_t i=0; i<m_neq; i++) {
667 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
668 weightedErrors.emplace_back(-abs(err), err, i);
673 N = std::min(N,
static_cast<int>(m_neq));
674 sort(weightedErrors.begin(), weightedErrors.end());
675 fmt::memory_buffer s;
676 for (
int i=0; i<N; i++) {
678 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
684 const string& cvodesMethod)
const
686 if (flag == CV_SUCCESS) {
688 }
else if (flag == CV_MEM_NULL) {
690 "CVODES integrator is not initialized");
692 const char* flagname = CVodeGetReturnFlagName(flag);
694 "{} 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< PreconditionerBase > 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.
void fmt_append(fmt::memory_buffer &b, Args... args)
Versions 6.2.0 and 6.2.1 of fmtlib do not include this define before they include windows....
static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
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.