9 #include "sundials/sundials_types.h"
10 #include "sundials/sundials_math.h"
12 #if CT_SUNDIALS_VERSION >= 30
13 #if CT_SUNDIALS_USE_LAPACK
14 #include "sunlinsol/sunlinsol_lapackdense.h"
15 #include "sunlinsol/sunlinsol_lapackband.h"
17 #include "sunlinsol/sunlinsol_dense.h"
18 #include "sunlinsol/sunlinsol_band.h"
20 #include "sunlinsol/sunlinsol_spgmr.h"
21 #include "ida/ida_direct.h"
22 #include "ida/ida_spils.h"
24 #include "ida/ida_dense.h"
25 #include "ida/ida_spgmr.h"
26 #include "ida/ida_band.h"
28 #include "nvector/nvector_serial.h"
32 #if CT_SUNDIALS_VERSION < 25
33 typedef int sd_size_t;
35 typedef long int sd_size_t;
78 static int ida_resid(realtype t, N_Vector y, N_Vector ydot, N_Vector r,
void* f_data)
85 int flag = f->
evalResidNJ(t, delta_t, NV_DATA_S(y), NV_DATA_S(ydot),
110 #if CT_SUNDIALS_VERSION >= 30
111 static int ida_jacobian(realtype t, realtype c_j, N_Vector y, N_Vector yp,
112 N_Vector r, SUNMatrix Jac,
void *f_data,
113 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
120 if (SUNMatGetID(Jac) == SUNMATRIX_DENSE) {
121 cols = SM_COLS_D(Jac);
122 }
else if (SUNMatGetID(Jac) == SUNMATRIX_BAND) {
123 cols = SM_COLS_B(Jac);
132 static int ida_jacobian(sd_size_t nrows, realtype t, realtype c_j, N_Vector y, N_Vector ydot, N_Vector r,
133 DlsMat Jac,
void* f_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
140 Jac->cols, NV_DATA_S(r));
176 m_maxErrTestFails(-1),
178 m_maxNonlinConvFails(-1),
185 IDA_Solver::~IDA_Solver()
191 N_VDestroy_Serial(
m_y);
194 N_VDestroy_Serial(
m_ydot);
197 N_VDestroy_Serial(m_abstol);
200 N_VDestroy_Serial(m_constraints);
206 return NV_Ith_S(
m_y,k);
209 const doublereal* IDA_Solver::solutionVector()
const
211 return NV_DATA_S(
m_y);
216 return NV_Ith_S(
m_ydot,k);
219 const doublereal* IDA_Solver::derivativeVector()
const
228 m_abstol = N_VNew_Serial(
m_neq);
230 for (
int i = 0; i <
m_neq; i++) {
231 NV_Ith_S(m_abstol, i) = abstol[i];
235 int flag = IDASVtolerances(
m_ida_mem, m_reltol, m_abstol);
236 if (flag != IDA_SUCCESS) {
238 "Memory allocation failed.");
249 int flag = IDASStolerances(
m_ida_mem, m_reltol, m_abstols);
250 if (flag != IDA_SUCCESS) {
252 "Memory allocation failed.");
257 void IDA_Solver::setLinearSolverType(
int solverType)
264 setLinearSolverType(0);
274 void IDA_Solver::setMaxOrder(
int n)
305 #if CT_SUNDIALS_VERSION >= 30
310 if (flag != IDA_SUCCESS) {
312 "IDADlsSetDenseJacFn failed.");
317 void IDA_Solver::setMaxErrTestFailures(
int maxErrTestFails)
332 void IDA_Solver::inclAlgebraicInErrorTest(
bool yesno)
348 N_VDestroy_Serial(
m_y);
351 N_VDestroy_Serial(
m_ydot);
354 N_VDestroy_Serial(m_id);
357 N_VDestroy_Serial(m_constraints);
362 m_constraints = N_VNew_Serial(
m_neq);
364 for (
int i=0; i<
m_neq; i++) {
365 NV_Ith_S(
m_y, i) = 0.0;
366 NV_Ith_S(
m_ydot, i) = 0.0;
367 NV_Ith_S(m_constraints, i) = 0.0;
381 if (m_itol == IDA_SV) {
383 if (flag != IDA_SUCCESS) {
384 if (flag == IDA_MEM_FAIL) {
386 "Memory allocation failed.");
387 }
else if (flag == IDA_ILL_INPUT) {
389 "Illegal value for IDAInit input argument.");
391 throw CanteraError(
"IDA_Solver::init",
"IDAInit failed.");
394 flag = IDASVtolerances(
m_ida_mem, m_reltol, m_abstol);
395 if (flag != IDA_SUCCESS) {
396 throw CanteraError(
"IDA_Solver::init",
"Memory allocation failed.");
400 if (flag != IDA_SUCCESS) {
401 if (flag == IDA_MEM_FAIL) {
403 "Memory allocation failed.");
404 }
else if (flag == IDA_ILL_INPUT) {
406 "Illegal value for IDAInit input argument.");
408 throw CanteraError(
"IDA_Solver::init",
"IDAInit failed.");
411 flag = IDASStolerances(
m_ida_mem, m_reltol, m_abstols);
412 if (flag != IDA_SUCCESS) {
413 throw CanteraError(
"IDA_Solver::init",
"Memory allocation failed.");
418 if (m_type == 1 || m_type == 0) {
421 #if CT_SUNDIALS_VERSION >= 30
422 SUNLinSolFree((SUNLinearSolver)
m_linsol);
427 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
429 #if CT_SUNDIALS_USE_LAPACK
440 throw CanteraError(
"IDA_Solver::init",
"IDADense failed");
442 }
else if (m_type == 2) {
444 long int nu = m_mupper;
445 long int nl = m_mlower;
446 #if CT_SUNDIALS_VERSION >= 30
447 SUNLinSolFree((SUNLinearSolver)
m_linsol);
449 #if CT_SUNDIALS_VERSION < 40
456 "Unable to create SUNBandMatrix of size {} with bandwidths "
457 "{} and {}", N, nu, nl);
459 #if CT_SUNDIALS_USE_LAPACK
471 "unsupported linear solver type");
475 #if CT_SUNDIALS_VERSION >= 30
480 if (flag != IDA_SUCCESS) {
482 "IDADlsSetDenseJacFn failed.");
488 flag = IDASetUserData(
m_ida_mem, m_fdata.get());
489 if (flag != IDA_SUCCESS) {
490 throw CanteraError(
"IDA_Solver::init",
"IDASetUserData failed.");
496 if (flag != IDA_SUCCESS) {
497 throw CanteraError(
"IDA_Solver::init",
"IDASetMaxOrd failed.");
502 if (flag != IDA_SUCCESS) {
503 throw CanteraError(
"IDA_Solver::init",
"IDASetMaxNumSteps failed.");
508 if (flag != IDA_SUCCESS) {
509 throw CanteraError(
"IDA_Solver::init",
"IDASetInitStep failed.");
514 if (flag != IDA_SUCCESS) {
515 throw CanteraError(
"IDA_Solver::init",
"IDASetStopTime failed.");
520 if (flag != IDA_SUCCESS) {
522 "IDASetMaxErrTestFails failed.");
527 if (flag != IDA_SUCCESS) {
529 "IDASetmaxNonlinIters failed.");
534 if (flag != IDA_SUCCESS) {
536 "IDASetMaxConvFails failed.");
541 if (flag != IDA_SUCCESS) {
542 throw CanteraError(
"IDA_Solver::init",
"IDASetSuppressAlg failed.");
549 doublereal tout1 = tout;
558 int flag = IDACalcIC(
m_ida_mem, IDA_Y_INIT, tout1);
559 if (flag != IDA_SUCCESS) {
560 throw CanteraError(
"IDA_Solver::correctInitial_Y_given_Yp",
561 "IDACalcIC failed: error = {}", flag);
565 if (flag != IDA_SUCCESS) {
566 throw CanteraError(
"IDA_Solver::correctInitial_Y_given_Yp",
567 "IDAGetSolution failed: error = {}", flag);
569 for (
int i = 0; i <
m_neq; i++) {
570 y[i] = NV_Ith_S(
m_y, i);
571 yp[i] = NV_Ith_S(
m_ydot, i);
577 int icopt = IDA_YA_YDP_INIT;
578 doublereal tout1 = tout;
587 int flag = IDACalcIC(
m_ida_mem, icopt, tout1);
588 if (flag != IDA_SUCCESS) {
589 throw CanteraError(
"IDA_Solver::correctInitial_YaYp_given_Yd",
590 "IDACalcIC failed: error = {}", flag);
594 if (flag != IDA_SUCCESS) {
595 throw CanteraError(
"IDA_Solver::correctInitial_YaYp_given_Yd",
596 "IDAGetSolution failed: error = {}", flag);
598 for (
int i = 0; i <
m_neq; i++) {
599 y[i] = NV_Ith_S(
m_y, i);
600 yp[i] = NV_Ith_S(
m_ydot, i);
606 double tretn = tout - 1000;
609 if (flag != IDA_SUCCESS) {
610 throw CanteraError(
"IDA_Solver::solve",
"IDA error encountered.");
612 while (tretn < tout) {
614 throw CanteraError(
"IDA_Solver::solve",
"tout <= tcurrent");
620 throw CanteraError(
"IDA_Solver::solve",
"IDA error encountered.");
621 }
else if (flag == IDA_TSTOP_RETURN) {
623 }
else if (flag == IDA_ROOT_RETURN) {
625 }
else if (flag == IDA_WARNING) {
626 throw CanteraError(
"IDA_Solver::solve",
"IDA Warning encountered.");
632 if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
633 throw CanteraError(
"IDA_Solver::solve",
"IDA error encountered.");
643 throw CanteraError(
"IDA_Solver::step",
"tout <= tcurrent");
649 throw CanteraError(
"IDA_Solver::step",
"IDA error encountered.");
650 }
else if (flag == IDA_TSTOP_RETURN) {
652 }
else if (flag == IDA_ROOT_RETURN) {
654 }
else if (flag == IDA_WARNING) {
655 throw CanteraError(
"IDA_Solver::step",
"IDA Warning encountered.");
664 long int lenrw, leniw;
666 case REAL_WORKSPACE_SIZE:
667 flag = IDAGetWorkSpace(
m_ida_mem, &lenrw, &leniw);
668 return doublereal(lenrw);
static int ida_jacobian(sd_size_t nrows, realtype t, realtype c_j, N_Vector y, N_Vector ydot, N_Vector r, DlsMat Jac, void *f_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Function called by by IDA to evaluate the Jacobian, given y and ydot.
static int ida_resid(realtype 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.
Header file for class IDA_Solver.
Base class for exceptions thrown by Cantera classes.
integer m_neq
Number of total equations in the system.
Wrapper for Sundials IDA solver.
virtual doublereal step(doublereal tout)
Take one internal step.
virtual void setStopTime(doublereal tstop)
Set the stop time.
virtual void correctInitial_Y_given_Yp(doublereal *y, doublereal *yp, doublereal tout)
Calculate consistent value of the starting solution given the starting solution derivatives.
virtual void setMaxNumSteps(int n)
Set the maximum number of time steps.
virtual void setBandedLinearSolver(int m_upper, int m_lower)
Set up the problem to use a band solver.
virtual doublereal solution(int k) const
the current value of solution component k.
N_Vector m_y
Current value of the solution vector.
virtual void setJacobianType(int formJac)
Set the form of the Jacobian.
virtual void correctInitial_YaYp_given_Yd(doublereal *y, doublereal *yp, doublereal tout)
Calculate consistent value of the algebraic constraints and derivatives at the start of the problem.
void * m_linsol
Sundials linear solver object.
int m_maxNonlinConvFails
Maximum number of nonlinear convergence failures.
virtual void setInitialStepSize(doublereal h0)
Set the initial step size.
virtual int solve(doublereal tout)
Step the system to a final value of the time.
int m_maxNonlinIters
Maximum number of nonlinear solver iterations at one solution.
doublereal m_told
Value of the previous time.
virtual void setDenseLinearSolver()
Set up the problem to use a dense linear direct solver.
doublereal m_deltat
Value of deltaT for the current step.
int m_maxErrTestFails
maximum number of error test failures
virtual void setTolerances(doublereal reltol, doublereal *abstol)
Set error tolerances.
virtual void setMaxNonlinIterations(int n)
Set the maximum number of nonlinear iterations on a timestep.
doublereal m_tcurrent
Value of the current time.
int m_maxord
maximum time step order of the method
virtual doublereal getOutputParameter(int flag) const
Get the value of a solver-specific output parameter.
virtual void init(doublereal t0)
initialize.
int m_formJac
Form of the Jacobian.
doublereal m_told_old
Value of the previous, previous time.
virtual void setMaxNonlinConvFailures(int n)
Set the maximum number of nonlinear solver convergence failures.
virtual double getCurrentStepFromIDA()
Get the current step size from IDA via a call.
virtual doublereal derivative(int k) const
the current value of the derivative of solution component k.
doublereal m_h0
Value of the initial time step.
void * m_linsol_matrix
matrix used by Sundials
N_Vector m_ydot
Current value of the derivative of the solution vector.
doublereal m_t0
Initial value of the time.
int m_maxsteps
Maximum number of time steps allowed.
doublereal m_tstop
maximum time
int m_setSuppressAlg
If true, the algebraic variables don't contribute to error tolerances.
void * m_ida_mem
Pointer to the IDA memory for the problem.
A simple class to hold an array of parameter values and a pointer to an instance of a subclass of Res...
int nparams() const
Return the number of parameters in the calculation.
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
virtual int evalJacobianDP(const doublereal t, const doublereal delta_t, doublereal cj, const doublereal *const y, const doublereal *const ydot, doublereal *const *jacobianColPts, doublereal *const resid)
Calculate an analytical Jacobian and the residual at the current time and values.
virtual int evalResidNJ(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot, doublereal *const resid, const ResidEval_Type_Enum evalType=Base_ResidEval, const int id_x=-1, const doublereal delta_x=0.0)
Evaluate the residual function.
virtual int getInitialConditions(const doublereal t0, doublereal *const y, doublereal *const ydot)
Fill in the initial conditions.
Namespace for the Cantera kernel.
Contains declarations for string manipulation functions within Cantera.