10 #if SUNDIALS_VERSION >= 24
11 #include "sundials/sundials_types.h"
12 #include "sundials/sundials_math.h"
14 #include "ida/ida_dense.h"
15 #include "ida/ida_spgmr.h"
16 #include "ida/ida_band.h"
17 #include "nvector/nvector_serial.h"
21 #if SUNDIALS_VERSION < 25
22 typedef int sd_size_t;
24 typedef long int sd_size_t;
68 static int ida_resid(realtype t, N_Vector y, N_Vector ydot, N_Vector r,
void* f_data)
70 double* ydata = NV_DATA_S(y);
71 double* ydotdata = NV_DATA_S(ydot);
72 double* rdata = NV_DATA_S(r);
75 Cantera::IDA_Solver* s = d->m_solver;
76 double delta_t = s->getCurrentStepFromIDA();
79 int flag = f->
evalResidNJ(t, delta_t, ydata, ydotdata, rdata);
101 static int ida_jacobian(sd_size_t nrows, realtype t, realtype c_j, N_Vector y, N_Vector ydot, N_Vector r,
102 DlsMat Jac,
void* f_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
104 doublereal* ydata = NV_DATA_S(y);
105 doublereal* ydotdata = NV_DATA_S(ydot);
106 doublereal* rdata = NV_DATA_S(r);
109 doublereal*
const* colPts = Jac->cols;
110 Cantera::IDA_Solver* s = d->m_solver;
111 double delta_t = s->getCurrentStepFromIDA();
113 f->
evalJacobianDP(t, delta_t, c_j, ydata, ydotdata, colPts, rdata);
121 IDA_Solver::IDA_Solver(ResidJacEval& f) :
147 m_maxErrTestFails(-1),
149 m_maxNonlinConvFails(-1),
157 IDA_Solver::~IDA_Solver()
163 N_VDestroy_Serial(m_y);
166 N_VDestroy_Serial(m_ydot);
169 N_VDestroy_Serial(m_abstol);
172 N_VDestroy_Serial(m_constraints);
177 doublereal IDA_Solver::solution(
int k)
const
179 return NV_Ith_S(m_y,k);
182 const doublereal* IDA_Solver::solutionVector()
const
184 return NV_DATA_S(m_y);
187 doublereal IDA_Solver::derivative(
int k)
const
189 return NV_Ith_S(m_ydot,k);
192 const doublereal* IDA_Solver::derivativeVector()
const
194 return NV_DATA_S(m_ydot);
197 void IDA_Solver::setTolerances(
double reltol,
double* abstol)
201 m_abstol = N_VNew_Serial(m_neq);
203 for (
int i = 0; i < m_neq; i++) {
204 NV_Ith_S(m_abstol, i) = abstol[i];
208 int flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
209 if (flag != IDA_SUCCESS) {
210 throw IDA_Err(
"Memory allocation failed.");
215 void IDA_Solver::setTolerances(doublereal reltol, doublereal abstol)
221 int flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
222 if (flag != IDA_SUCCESS) {
223 throw IDA_Err(
"Memory allocation failed.");
228 void IDA_Solver::setLinearSolverType(
int solverType)
233 void IDA_Solver::setDenseLinearSolver()
235 setLinearSolverType(0);
238 void IDA_Solver::setBandedLinearSolver(
int m_upper,
int m_lower)
245 void IDA_Solver::setMaxOrder(
int n)
250 void IDA_Solver::setMaxNumSteps(
int n)
255 void IDA_Solver::setInitialStepSize(doublereal h0)
260 void IDA_Solver::setStopTime(doublereal tstop)
265 doublereal IDA_Solver::getCurrentStepFromIDA()
268 IDAGetCurrentStep(m_ida_mem, &hcur);
272 void IDA_Solver::setJacobianType(
int formJac)
276 if (m_formJac == 1) {
277 int flag = IDADlsSetDenseJacFn(m_ida_mem,
ida_jacobian);
278 if (flag != IDA_SUCCESS) {
279 throw IDA_Err(
"IDADlsSetDenseJacFn failed.");
285 void IDA_Solver::setMaxErrTestFailures(
int maxErrTestFails)
287 m_maxErrTestFails = maxErrTestFails;
290 void IDA_Solver::setMaxNonlinIterations(
int n)
292 m_maxNonlinIters = n;
295 void IDA_Solver::setMaxNonlinConvFailures(
int n)
297 m_maxNonlinConvFails = n;
300 void IDA_Solver::inclAlgebraicInErrorTest(
bool yesno)
303 m_setSuppressAlg = 0;
305 m_setSuppressAlg = 1;
309 void IDA_Solver::init(doublereal t0)
316 N_VDestroy_Serial(m_y);
319 N_VDestroy_Serial(m_ydot);
322 N_VDestroy_Serial(m_id);
325 N_VDestroy_Serial(m_constraints);
328 m_y = N_VNew_Serial(m_neq);
329 m_ydot = N_VNew_Serial(m_neq);
330 m_constraints = N_VNew_Serial(m_neq);
332 for (
int i=0; i<m_neq; i++) {
333 NV_Ith_S(m_y, i) = 0.0;
334 NV_Ith_S(m_ydot, i) = 0.0;
335 NV_Ith_S(m_constraints, i) = 0.0;
339 m_resid.getInitialConditions(m_t0, NV_DATA_S(m_y), NV_DATA_S(m_ydot));
346 m_ida_mem = IDACreate();
352 if (m_itol == IDA_SV) {
353 #if SUNDIALS_VERSION <= 23
355 flag = IDAMalloc(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot,
356 m_itol, m_reltol, m_abstol);
357 if (flag != IDA_SUCCESS) {
358 if (flag == IDA_MEM_FAIL) {
359 throw IDA_Err(
"Memory allocation failed.");
360 }
else if (flag == IDA_ILL_INPUT) {
361 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
363 throw IDA_Err(
"IDAMalloc failed.");
367 #elif SUNDIALS_VERSION >= 24
368 flag = IDAInit(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot);
369 if (flag != IDA_SUCCESS) {
370 if (flag == IDA_MEM_FAIL) {
371 throw IDA_Err(
"Memory allocation failed.");
372 }
else if (flag == IDA_ILL_INPUT) {
373 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
375 throw IDA_Err(
"IDAMalloc failed.");
378 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
379 if (flag != IDA_SUCCESS) {
380 throw IDA_Err(
"Memory allocation failed.");
384 #if SUNDIALS_VERSION <= 23
386 flag = IDAMalloc(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot,
387 m_itol, m_reltol, &m_abstols);
388 if (flag != IDA_SUCCESS) {
389 if (flag == IDA_MEM_FAIL) {
390 throw IDA_Err(
"Memory allocation failed.");
391 }
else if (flag == IDA_ILL_INPUT) {
392 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
394 throw IDA_Err(
"IDAMalloc failed.");
398 #elif SUNDIALS_VERSION >= 24
399 flag = IDAInit(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot);
400 if (flag != IDA_SUCCESS) {
401 if (flag == IDA_MEM_FAIL) {
402 throw IDA_Err(
"Memory allocation failed.");
403 }
else if (flag == IDA_ILL_INPUT) {
404 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
406 throw IDA_Err(
"IDAMalloc failed.");
409 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
410 if (flag != IDA_SUCCESS) {
411 throw IDA_Err(
"Memory allocation failed.");
420 if (m_type == 1 || m_type == 0) {
422 flag = IDADense(m_ida_mem, N);
424 throw IDA_Err(
"IDADense failed");
426 }
else if (m_type == 2) {
428 long int nu = m_mupper;
429 long int nl = m_mlower;
430 IDABand(m_ida_mem, N, nu, nl);
432 throw IDA_Err(
"unsupported linear solver type");
435 if (m_formJac == 1) {
437 if (flag != IDA_SUCCESS) {
438 throw IDA_Err(
"IDADlsSetDenseJacFn failed.");
443 m_fdata =
new ResidData(&m_resid,
this, m_resid.nparams());
444 #if SUNDIALS_VERSION <= 23
445 flag = IDASetRdata(m_ida_mem, (
void*)m_fdata);
446 if (flag != IDA_SUCCESS) {
447 throw IDA_Err(
"IDASetRdata failed.");
449 #elif SUNDIALS_VERSION >= 24
450 flag = IDASetUserData(m_ida_mem, (
void*)m_fdata);
451 if (flag != IDA_SUCCESS) {
452 throw IDA_Err(
"IDASetUserData failed.");
458 flag = IDASetMaxOrd(m_ida_mem, m_maxord);
459 if (flag != IDA_SUCCESS) {
460 throw IDA_Err(
"IDASetMaxOrd failed.");
463 if (m_maxsteps > 0) {
464 flag = IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
465 if (flag != IDA_SUCCESS) {
466 throw IDA_Err(
"IDASetMaxNumSteps failed.");
470 flag = IDASetInitStep(m_ida_mem, m_h0);
471 if (flag != IDA_SUCCESS) {
472 throw IDA_Err(
"IDASetInitStep failed.");
476 flag = IDASetStopTime(m_ida_mem, m_tstop);
477 if (flag != IDA_SUCCESS) {
478 throw IDA_Err(
"IDASetStopTime failed.");
481 if (m_maxErrTestFails >= 0) {
482 flag = IDASetMaxErrTestFails(m_ida_mem, m_maxErrTestFails);
483 if (flag != IDA_SUCCESS) {
484 throw IDA_Err(
"IDASetMaxErrTestFails failed.");
487 if (m_maxNonlinIters > 0) {
488 flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
489 if (flag != IDA_SUCCESS) {
490 throw IDA_Err(
"IDASetmaxNonlinIters failed.");
493 if (m_maxNonlinConvFails >= 0) {
494 flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
495 if (flag != IDA_SUCCESS) {
496 throw IDA_Err(
"IDASetMaxConvFails failed.");
499 if (m_setSuppressAlg != 0) {
500 flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
501 if (flag != IDA_SUCCESS) {
502 throw IDA_Err(
"IDASetSuppressAlg failed.");
507 void IDA_Solver::correctInitial_Y_given_Yp(doublereal* y, doublereal* yp, doublereal tout)
509 int icopt = IDA_Y_INIT;
510 doublereal tout1 = tout;
519 int flag = IDACalcIC(m_ida_mem, icopt, tout1);
520 if (flag != IDA_SUCCESS) {
521 throw IDA_Err(
"IDACalcIC failed: error = " +
int2str(flag));
525 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
526 if (flag != IDA_SUCCESS) {
527 throw IDA_Err(
"IDAGetSolution failed: error = " +
int2str(flag));
529 doublereal* yy = NV_DATA_S(m_y);
530 doublereal* yyp = NV_DATA_S(m_ydot);
532 for (
int i = 0; i < m_neq; i++) {
538 void IDA_Solver::correctInitial_YaYp_given_Yd(doublereal* y, doublereal* yp, doublereal tout)
540 int icopt = IDA_YA_YDP_INIT;
541 doublereal tout1 = tout;
550 int flag = IDACalcIC(m_ida_mem, icopt, tout1);
551 if (flag != IDA_SUCCESS) {
552 throw IDA_Err(
"IDACalcIC failed: error = " +
int2str(flag));
556 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
557 if (flag != IDA_SUCCESS) {
558 throw IDA_Err(
"IDAGetSolution failed: error = " +
int2str(flag));
560 doublereal* yy = NV_DATA_S(m_y);
561 doublereal* yyp = NV_DATA_S(m_ydot);
563 for (
int i = 0; i < m_neq; i++) {
573 flag = IDASetStopTime(m_ida_mem, tout);
574 if (flag != IDA_SUCCESS) {
575 throw IDA_Err(
" IDA error encountered.");
578 if (tout <= m_tcurrent) {
579 throw IDA_Err(
" tout <= tcurrent");
583 flag = IDASolve(m_ida_mem, tout, &tretn, m_y, m_ydot, IDA_ONE_STEP);
585 throw IDA_Err(
" IDA error encountered.");
586 }
else if (flag == IDA_TSTOP_RETURN) {
588 }
else if (flag == IDA_ROOT_RETURN) {
590 }
else if (flag == IDA_WARNING) {
591 throw IDA_Err(
" IDA Warning encountered.");
594 m_deltat = m_tcurrent - m_told;
595 }
while (tretn < tout);
597 if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
598 throw IDA_Err(
" IDA error encountered.");
603 double IDA_Solver::step(
double tout)
607 if (tout <= m_tcurrent) {
608 throw IDA_Err(
" tout <= tcurrent");
612 flag = IDASolve(m_ida_mem, tout, &t, m_y, m_ydot, IDA_ONE_STEP);
614 throw IDA_Err(
" IDA error encountered.");
615 }
else if (flag == IDA_TSTOP_RETURN) {
617 }
else if (flag == IDA_ROOT_RETURN) {
619 }
else if (flag == IDA_WARNING) {
620 throw IDA_Err(
" IDA Warning encountered.");
623 m_deltat = m_tcurrent - m_told;
627 doublereal IDA_Solver::getOutputParameter(
int flag)
const
629 long int lenrw, leniw;
631 case REAL_WORKSPACE_SIZE:
632 flag = IDAGetWorkSpace(m_ida_mem, &lenrw, &leniw);
633 return doublereal(lenrw);
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. ...
A simple class to hold an array of parameter values and a pointer to an instance of a subclass of Res...
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
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.
int solve(DenseMatrix &A, double *b)
Solve Ax = b. Array b is overwritten on exit with x.
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.
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
Contains declarations for string manipulation functions within Cantera.
Header file for class IDA_Solver.
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.