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();
112 f->
evalJacobianDP(t, delta_t, c_j, ydata, ydotdata, colPts, rdata);
120 IDA_Solver::IDA_Solver(ResidJacEval& f) :
146 m_maxErrTestFails(-1),
148 m_maxNonlinConvFails(-1),
156 IDA_Solver::~IDA_Solver()
162 N_VDestroy_Serial(m_y);
165 N_VDestroy_Serial(m_ydot);
168 N_VDestroy_Serial(m_abstol);
171 N_VDestroy_Serial(m_constraints);
176 doublereal IDA_Solver::solution(
int k)
const
178 return NV_Ith_S(m_y,k);
181 const doublereal* IDA_Solver::solutionVector()
const
183 return NV_DATA_S(m_y);
186 doublereal IDA_Solver::derivative(
int k)
const
188 return NV_Ith_S(m_ydot,k);
191 const doublereal* IDA_Solver::derivativeVector()
const
193 return NV_DATA_S(m_ydot);
196 void IDA_Solver::setTolerances(
double reltol,
double* abstol)
200 m_abstol = N_VNew_Serial(m_neq);
202 for (
int i = 0; i < m_neq; i++) {
203 NV_Ith_S(m_abstol, i) = abstol[i];
207 int flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
208 if (flag != IDA_SUCCESS) {
209 throw IDA_Err(
"Memory allocation failed.");
214 void IDA_Solver::setTolerances(doublereal reltol, doublereal abstol)
220 int flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
221 if (flag != IDA_SUCCESS) {
222 throw IDA_Err(
"Memory allocation failed.");
227 void IDA_Solver::setLinearSolverType(
int solverType)
232 void IDA_Solver::setDenseLinearSolver()
234 setLinearSolverType(0);
237 void IDA_Solver::setBandedLinearSolver(
int m_upper,
int m_lower)
244 void IDA_Solver::setMaxOrder(
int n)
249 void IDA_Solver::setMaxNumSteps(
int n)
254 void IDA_Solver::setInitialStepSize(doublereal h0)
259 void IDA_Solver::setStopTime(doublereal tstop)
264 doublereal IDA_Solver::getCurrentStepFromIDA()
267 IDAGetCurrentStep(m_ida_mem, &hcur);
271 void IDA_Solver::setJacobianType(
int formJac)
275 if (m_formJac == 1) {
276 int flag = IDADlsSetDenseJacFn(m_ida_mem,
ida_jacobian);
277 if (flag != IDA_SUCCESS) {
278 throw IDA_Err(
"IDADlsSetDenseJacFn failed.");
284 void IDA_Solver::setMaxErrTestFailures(
int maxErrTestFails)
286 m_maxErrTestFails = maxErrTestFails;
289 void IDA_Solver::setMaxNonlinIterations(
int n)
291 m_maxNonlinIters = n;
294 void IDA_Solver::setMaxNonlinConvFailures(
int n)
296 m_maxNonlinConvFails = n;
299 void IDA_Solver::inclAlgebraicInErrorTest(
bool yesno)
302 m_setSuppressAlg = 0;
304 m_setSuppressAlg = 1;
308 void IDA_Solver::init(doublereal t0)
315 N_VDestroy_Serial(m_y);
318 N_VDestroy_Serial(m_ydot);
321 N_VDestroy_Serial(m_id);
324 N_VDestroy_Serial(m_constraints);
327 m_y = N_VNew_Serial(m_neq);
328 m_ydot = N_VNew_Serial(m_neq);
329 m_constraints = N_VNew_Serial(m_neq);
331 for (
int i=0; i<m_neq; i++) {
332 NV_Ith_S(m_y, i) = 0.0;
333 NV_Ith_S(m_ydot, i) = 0.0;
334 NV_Ith_S(m_constraints, i) = 0.0;
338 m_resid.getInitialConditions(m_t0, NV_DATA_S(m_y), NV_DATA_S(m_ydot));
345 m_ida_mem = IDACreate();
351 if (m_itol == IDA_SV) {
352 flag = IDAInit(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot);
353 if (flag != IDA_SUCCESS) {
354 if (flag == IDA_MEM_FAIL) {
355 throw IDA_Err(
"Memory allocation failed.");
356 }
else if (flag == IDA_ILL_INPUT) {
357 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
359 throw IDA_Err(
"IDAMalloc failed.");
362 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
363 if (flag != IDA_SUCCESS) {
364 throw IDA_Err(
"Memory allocation failed.");
367 flag = IDAInit(m_ida_mem,
ida_resid, m_t0, m_y, m_ydot);
368 if (flag != IDA_SUCCESS) {
369 if (flag == IDA_MEM_FAIL) {
370 throw IDA_Err(
"Memory allocation failed.");
371 }
else if (flag == IDA_ILL_INPUT) {
372 throw IDA_Err(
"Illegal value for IDAMalloc input argument.");
374 throw IDA_Err(
"IDAMalloc failed.");
377 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
378 if (flag != IDA_SUCCESS) {
379 throw IDA_Err(
"Memory allocation failed.");
387 if (m_type == 1 || m_type == 0) {
389 flag = IDADense(m_ida_mem, N);
391 throw IDA_Err(
"IDADense failed");
393 }
else if (m_type == 2) {
395 long int nu = m_mupper;
396 long int nl = m_mlower;
397 IDABand(m_ida_mem, N, nu, nl);
399 throw IDA_Err(
"unsupported linear solver type");
402 if (m_formJac == 1) {
404 if (flag != IDA_SUCCESS) {
405 throw IDA_Err(
"IDADlsSetDenseJacFn failed.");
410 m_fdata =
new ResidData(&m_resid,
this, m_resid.nparams());
411 flag = IDASetUserData(m_ida_mem, (
void*)m_fdata);
412 if (flag != IDA_SUCCESS) {
413 throw IDA_Err(
"IDASetUserData failed.");
418 flag = IDASetMaxOrd(m_ida_mem, m_maxord);
419 if (flag != IDA_SUCCESS) {
420 throw IDA_Err(
"IDASetMaxOrd failed.");
423 if (m_maxsteps > 0) {
424 flag = IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
425 if (flag != IDA_SUCCESS) {
426 throw IDA_Err(
"IDASetMaxNumSteps failed.");
430 flag = IDASetInitStep(m_ida_mem, m_h0);
431 if (flag != IDA_SUCCESS) {
432 throw IDA_Err(
"IDASetInitStep failed.");
436 flag = IDASetStopTime(m_ida_mem, m_tstop);
437 if (flag != IDA_SUCCESS) {
438 throw IDA_Err(
"IDASetStopTime failed.");
441 if (m_maxErrTestFails >= 0) {
442 flag = IDASetMaxErrTestFails(m_ida_mem, m_maxErrTestFails);
443 if (flag != IDA_SUCCESS) {
444 throw IDA_Err(
"IDASetMaxErrTestFails failed.");
447 if (m_maxNonlinIters > 0) {
448 flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
449 if (flag != IDA_SUCCESS) {
450 throw IDA_Err(
"IDASetmaxNonlinIters failed.");
453 if (m_maxNonlinConvFails >= 0) {
454 flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
455 if (flag != IDA_SUCCESS) {
456 throw IDA_Err(
"IDASetMaxConvFails failed.");
459 if (m_setSuppressAlg != 0) {
460 flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
461 if (flag != IDA_SUCCESS) {
462 throw IDA_Err(
"IDASetSuppressAlg failed.");
467 void IDA_Solver::correctInitial_Y_given_Yp(doublereal* y, doublereal* yp, doublereal tout)
469 int icopt = IDA_Y_INIT;
470 doublereal tout1 = tout;
479 int flag = IDACalcIC(m_ida_mem, icopt, tout1);
480 if (flag != IDA_SUCCESS) {
481 throw IDA_Err(
"IDACalcIC failed: error = " +
int2str(flag));
485 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
486 if (flag != IDA_SUCCESS) {
487 throw IDA_Err(
"IDAGetSolution failed: error = " +
int2str(flag));
489 doublereal* yy = NV_DATA_S(m_y);
490 doublereal* yyp = NV_DATA_S(m_ydot);
492 for (
int i = 0; i < m_neq; i++) {
498 void IDA_Solver::correctInitial_YaYp_given_Yd(doublereal* y, doublereal* yp, doublereal tout)
500 int icopt = IDA_YA_YDP_INIT;
501 doublereal tout1 = tout;
510 int flag = IDACalcIC(m_ida_mem, icopt, tout1);
511 if (flag != IDA_SUCCESS) {
512 throw IDA_Err(
"IDACalcIC failed: error = " +
int2str(flag));
516 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
517 if (flag != IDA_SUCCESS) {
518 throw IDA_Err(
"IDAGetSolution failed: error = " +
int2str(flag));
520 doublereal* yy = NV_DATA_S(m_y);
521 doublereal* yyp = NV_DATA_S(m_ydot);
523 for (
int i = 0; i < m_neq; i++) {
533 flag = IDASetStopTime(m_ida_mem, tout);
534 if (flag != IDA_SUCCESS) {
535 throw IDA_Err(
" IDA error encountered.");
538 if (tout <= m_tcurrent) {
539 throw IDA_Err(
" tout <= tcurrent");
543 flag = IDASolve(m_ida_mem, tout, &tretn, m_y, m_ydot, IDA_ONE_STEP);
545 throw IDA_Err(
" IDA error encountered.");
546 }
else if (flag == IDA_TSTOP_RETURN) {
548 }
else if (flag == IDA_ROOT_RETURN) {
550 }
else if (flag == IDA_WARNING) {
551 throw IDA_Err(
" IDA Warning encountered.");
554 m_deltat = m_tcurrent - m_told;
555 }
while (tretn < tout);
557 if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
558 throw IDA_Err(
" IDA error encountered.");
563 double IDA_Solver::step(
double tout)
567 if (tout <= m_tcurrent) {
568 throw IDA_Err(
" tout <= tcurrent");
572 flag = IDASolve(m_ida_mem, tout, &t, m_y, m_ydot, IDA_ONE_STEP);
574 throw IDA_Err(
" IDA error encountered.");
575 }
else if (flag == IDA_TSTOP_RETURN) {
577 }
else if (flag == IDA_ROOT_RETURN) {
579 }
else if (flag == IDA_WARNING) {
580 throw IDA_Err(
" IDA Warning encountered.");
583 m_deltat = m_tcurrent - m_told;
587 doublereal IDA_Solver::getOutputParameter(
int flag)
const
589 long int lenrw, leniw;
591 case REAL_WORKSPACE_SIZE:
592 flag = IDAGetWorkSpace(m_ida_mem, &lenrw, &leniw);
593 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, size_t nrhs, size_t ldb)
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.