9#include "cantera/numerics/sundials_headers.h"
17 return N_VNew_Serial(
static_cast<sd_size_t
>(N), context.get());
36static int ida_rhs(sunrealtype t, N_Vector y, N_Vector ydot, N_Vector r,
void* f_data)
42#if SUNDIALS_VERSION_MAJOR >= 7
47 static void sundials_err(
int line,
const char *func,
const char *file,
48 const char *msg, SUNErrCode err_code,
49 void *err_user_data, SUNContext sunctx)
51 IdasIntegrator* integrator = (IdasIntegrator*) err_user_data;
52 integrator->m_error_message = fmt::format(
"{}: {}\n", func, msg);
58 static void ida_err(
int error_code,
const char* module,
59 const char* function,
char* msg,
void* eh_data)
74IdasIntegrator::~IdasIntegrator()
80 N_VDestroy_Serial(
m_y);
89 N_VDestroy_Serial(m_constraints);
95 return NV_Ith_S(
m_y, k);
105 size_t n = abstol.size();
114 for (
size_t i=0; i<n; i++) {
143 checkError(flag,
"setMaxOrder",
"IDASetMaxOrd");
152 int flag = IDASetMaxStep(
m_ida_mem, hmax);
153 checkError(flag,
"setMaxStepSize",
"IDASetMaxStep");
184 long int steps = 0, stepSolveFails = 0, resEvals = 0, errTestFails = 0,
185 jacEvals = 0, linSetup = 0, linResEvals = 0, linIters = 0,
186 linConvFails = 0, precEvals = 0, precSolves = 0, jtSetupEvals = 0,
187 jTimesEvals = 0, nonlinIters = 0, nonlinConvFails = 0;
190 int flag = IDAGetNumSteps(
m_ida_mem, &steps);
191 checkError(flag,
"solverStats",
"IDAGetNumSteps");
195 IDAGetNumStepSolveFails(
m_ida_mem, &stepSolveFails);
197 IDAGetNumNonlinSolvIters(
m_ida_mem, &nonlinIters);
198 IDAGetNumNonlinSolvConvFails(
m_ida_mem, &nonlinConvFails);
199 IDAGetNumErrTestFails(
m_ida_mem, &errTestFails);
203 IDAGetNumLinResEvals(
m_ida_mem, &linResEvals);
204 IDAGetNumLinSolvSetups(
m_ida_mem, &linSetup);
206 IDAGetNumLinConvFails(
m_ida_mem, &linConvFails);
207 IDAGetNumPrecEvals(
m_ida_mem, &precEvals);
208 IDAGetNumPrecSolves(
m_ida_mem, &precSolves);
209 IDAGetNumJTSetupEvals(
m_ida_mem, &jtSetupEvals);
210 IDAGetNumJtimesEvals(
m_ida_mem, &jTimesEvals);
212 stats[
"steps"] = steps;
213 stats[
"step_solve_fails"] = stepSolveFails;
214 stats[
"res_evals"] = resEvals;
215 stats[
"rhs_evals"] = resEvals;
216 stats[
"nonlinear_iters"] = nonlinIters;
217 stats[
"nonlinear_conv_fails"] = nonlinConvFails;
218 stats[
"err_test_fails"] = errTestFails;
221 stats[
"jac_evals"] = jacEvals;
222 stats[
"lin_solve_setups"] = linSetup;
223 stats[
"lin_rhs_evals"] = linResEvals;
224 stats[
"lin_iters"] = linIters;
225 stats[
"lin_conv_fails"] = linConvFails;
226 stats[
"prec_evals"] = precEvals;
227 stats[
"prec_solves"] = precSolves;
228 stats[
"jt_vec_setup_evals"] = jtSetupEvals;
229 stats[
"jt_vec_prod_evals"] = jTimesEvals;
233void IdasIntegrator::setMaxNonlinIterations(
int n)
238 checkError(flag,
"setMaxNonlinIterations",
"IDASetMaxNonlinIters");
242void IdasIntegrator::setMaxNonlinConvFailures(
int n)
247 checkError(flag,
"setMaxNonlinConvFailures",
"IDASetMaxConvFails");
251void IdasIntegrator::includeAlgebraicInErrorTest(
bool yesno)
256 checkError(flag,
"inclAlgebraicInErrorTest",
"IDASetSuppressAlg");
270 N_VDestroy_Serial(
m_y);
276 N_VDestroy_Serial(
m_ydot);
284 "not enough absolute tolerance values specified.");
288 N_VDestroy_Serial(m_constraints);
304 throw CanteraError(
"IdasIntegrator::initialize",
"IDACreate failed.");
308 if (flag != IDA_SUCCESS) {
309 if (flag == IDA_MEM_FAIL) {
311 "Memory allocation failed.");
312 }
else if (flag == IDA_ILL_INPUT) {
314 "Illegal value for IDAInit input argument.");
316 throw CanteraError(
"IdasIntegrator::initialize",
"IDAInit failed.");
320 #if SUNDIALS_VERSION_MAJOR >= 7
321 flag = SUNContext_PushErrHandler(
m_sundials_ctx.get(), &sundials_err,
this);
327 flag = IDASetId(
m_ida_mem, m_constraints);
332 checkError(flag,
"initialize",
"IDASVtolerances");
335 checkError(flag,
"initialize",
"IDASStolerances");
339 checkError(flag,
"initialize",
"IDASetUserData");
342 throw CanteraError(
"IdasIntegrator::initialize",
"Sensitivity analysis "
343 "for DAE systems is not fully implemented");
347 checkError(flag,
"initialize",
"IDASetSensParams");
352void IdasIntegrator::reinitialize(
double t0,
FuncEval& func)
362 checkError(result,
"reinitialize",
"IDAReInit");
369 sd_size_t N =
static_cast<sd_size_t
>(
m_neq);
370 SUNLinSolFree((SUNLinearSolver)
m_linsol);
373 #if CT_SUNDIALS_USE_LAPACK
382 }
else if (
m_type ==
"GMRES") {
383 SUNLinSolFree((SUNLinearSolver)
m_linsol);
388 "unsupported linear solver flag '{}'",
m_type);
397 checkError(flag,
"applyOptions",
"IDASetMaxOrd");
407 checkError(flag,
"applyOptions",
"IDASetMaxNonlinIters");
411 checkError(flag,
"applyOptions",
"IDASetMaxConvFails");
415 checkError(flag,
"applyOptions",
"IDASetSuppressAlg");
419void IdasIntegrator::sensInit(
double t0,
FuncEval& func)
425 m_yS = N_VCloneVectorArray(
static_cast<int>(
m_np), y);
426 for (
size_t n = 0; n <
m_np; n++) {
427 N_VConst(0.0,
m_yS[n]);
429 N_VDestroy_Serial(y);
431 m_ySdot = N_VCloneVectorArray(
static_cast<int>(
m_np), ydot);
432 for (
size_t n = 0; n <
m_np; n++) {
436 int flag = IDASensInit(
m_ida_mem,
static_cast<sd_size_t
>(
m_np),
440 vector<double> atol(
m_np);
441 for (
size_t n = 0; n <
m_np; n++) {
447 checkError(flag,
"sensInit",
"IDASensSStolerances");
454 }
else if (tout <
m_time) {
456 "Cannot integrate backwards in time.\n"
457 "Requested time {} < current time {}",
464 "Maximum number of timesteps ({}) taken without reaching output "
465 "time ({}).\nCurrent integrator time: {}",
469 if (flag != IDA_SUCCESS) {
471 if (!f_errs.empty()) {
472 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
475 "IDA error encountered. Error code: {}\n{}\n"
477 "Components with largest weighted error estimates:\n{}",
490 if (flag != IDA_SUCCESS) {
492 if (!f_errs.empty()) {
493 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
496 "IDA error encountered. Error code: {}\n{}\n"
498 "Components with largest weighted error estimates:\n{}",
506double IdasIntegrator::sensitivity(
size_t k,
size_t p)
514 checkError(flag,
"sensitivity",
"IDAGetSens");
519 throw CanteraError(
"IdasIntegrator::sensitivity",
520 "sensitivity: k out of range ({})", k);
523 throw CanteraError(
"IdasIntegrator::sensitivity",
524 "sensitivity: p out of range ({})", p);
526 return NV_Ith_S(
m_yS[p],k);
536 vector<tuple<double, double, size_t>> weightedErrors;
537 for (
size_t i=0; i<
m_neq; i++) {
538 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
539 weightedErrors.emplace_back(-abs(err), err, i);
544 N = std::min(N,
static_cast<int>(
m_neq));
545 sort(weightedErrors.begin(), weightedErrors.end());
546 fmt::memory_buffer s;
547 for (
int i=0; i<N; i++) {
548 fmt_append(s,
"{}: {}\n", get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
554 const string& idaMethod)
const
556 if (flag == IDA_SUCCESS) {
558 }
else if (flag == IDA_MEM_NULL) {
560 "IDAS integrator is not initialized");
562 const char* flagname = IDAGetReturnFlagName(flag);
564 "{} returned error code {} ({}):\n{}",
573 throw CanteraError(
"IdasIntegrator::setMethod",
"unknown method");
Header file for class IdasIntegrator.
A map of string keys to values whose type can vary at runtime.
Base class for exceptions thrown by Cantera classes.
Virtual base class for ODE/DAE right-hand-side function evaluators.
virtual void getStateDae(span< double > y, span< double > ydot)
Fill in the vectors y and ydot with the current state of the system.
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
void clearErrors()
Clear any previously-stored suppressed errors.
virtual size_t neq() const =0
Number of equations.
int evalDaeNoThrow(double t, span< const double > y, span< const double > ydot, span< double > residual)
Evaluate the right-hand side using return code to indicate status.
virtual void getConstraints(span< double > constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
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...
Wrapper for Sundials IDAS solver.
span< double > solution() override
The current value of the solution of the system of equations.
double step(double tout) override
Integrate the system of equations.
void setMaxStepSize(double hmax) override
Set the maximum step size.
double m_t0
The start time for the integrator.
N_Vector m_y
The current system state.
void * m_linsol
Sundials linear solver object.
int m_maxNonlinConvFails
Maximum number of nonlinear convergence failures.
double m_init_step
Initial IDA step size.
void checkError(long flag, const string &ctMethod, const string &idaMethod) const
Check whether an IDAS method indicated an error.
double m_abstolsens
Scalar absolute tolerance for sensitivities.
size_t m_nabs
! Number of variables for which absolute tolerances were provided
void setMaxOrder(int n) override
Set the maximum integration order that will be used.
double m_time
The current integrator time, corresponding to m_y.
int m_maxNonlinIters
Maximum number of nonlinear solver iterations at one solution.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS and m_ySdot have been updated for the current inte...
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.
FuncEval * m_func
Object implementing the DAE residual function .
double m_abstols
Scalar absolute tolerance.
double m_hmax
Maximum time step size.
int m_itol
Flag indicating whether scalar (IDA_SS) or vector (IDA_SV) absolute tolerances are being used.
int m_maxErrTestFails
Maximum number of error test failures in attempting one step.
int m_maxord
Maximum order allowed for the BDF method.
void setSensitivityTolerances(double reltol, double abstol) override
Set the sensitivity error tolerances.
void applyOptions()
Applies user-specified options to the underlying IDAS solver.
void integrate(double tout) override
Integrate the system of equations.
void setLinearSolverType(const string &linearSolverType) override
Set the linear solver type.
double m_reltol
Relative tolerance for all state variables.
string m_error_message
Error message information provide by IDAS.
double m_reltolsens
Scalar relative tolerance for sensitivities.
N_Vector * m_yS
Sensitivities of y, size m_np by m_neq.
void setTolerances(double reltol, span< const double > abstol) override
Set error tolerances.
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.
size_t m_neq
Number of equations/variables in the system.
AnyMap solverStats() const override
Get solver stats from integrator.
string m_type
The linear solver type.
size_t m_np
Number of sensitivity parameters.
void * m_linsol_matrix
matrix used by Sundials
N_Vector m_ydot
The time derivatives of the system state.
N_Vector * m_ySdot
Sensitivities of ydot, size m_np by m_neq.
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.
int m_maxsteps
Maximum number of internal steps taken in a call to integrate()
string getErrorInfo(int N)
Returns a string listing the weighted error estimates associated with each solution component.
IdasIntegrator()
Constructor.
void setMethod(MethodType t) override
Set the solution method.
void * m_ida_mem
Pointer to the IDA memory for the problem.
bool m_setSuppressAlg
If true, the algebraic variables don't contribute to error tolerances.
N_Vector m_abstol
Absolute tolerances for each state variable.
virtual string linearSolverType() const
Return the integrator problem type.
virtual int lastOrder() const
Order used during the last solution step.
A wrapper for managing a SUNContext object.
Namespace for the Cantera kernel.
static int ida_rhs(sunrealtype t, N_Vector y, N_Vector ydot, N_Vector r, void *f_data)
Function called by IDA to evaluate the residual, given y and ydot.
static void ida_err(int error_code, const char *module, const char *function, char *msg, void *eh_data)
Function called by IDA 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.
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
Contains declarations for string manipulation functions within Cantera.