Cantera 2.6.0
IDA_Solver.cpp
Go to the documentation of this file.
1//! @file IDA_Solver.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
8
9#include "sundials/sundials_types.h"
10#include "sundials/sundials_math.h"
11#include "ida/ida.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"
16 #else
17 #include "sunlinsol/sunlinsol_dense.h"
18 #include "sunlinsol/sunlinsol_band.h"
19 #endif
20 #include "sunlinsol/sunlinsol_spgmr.h"
21 #include "ida/ida_direct.h"
22 #include "ida/ida_spils.h"
23#else
24 #include "ida/ida_dense.h"
25 #include "ida/ida_spgmr.h"
26 #include "ida/ida_band.h"
27#endif
28#include "nvector/nvector_serial.h"
29
30using namespace std;
31
32#if CT_SUNDIALS_VERSION < 25
33typedef int sd_size_t;
34#else
35typedef long int sd_size_t;
36#endif
37
38namespace Cantera
39{
40
41//! A simple class to hold an array of parameter values and a pointer to an
42//! instance of a subclass of ResidEval.
44{
45public:
46 ResidData(ResidJacEval* f, IDA_Solver* s, int npar = 0) {
47 m_func = f;
48 m_solver = s;
49 }
50
51 virtual ~ResidData() {
52 }
53
54 ResidJacEval* m_func;
55 IDA_Solver* m_solver;
56};
57}
58
59extern "C" {
60 //! Function called by IDA to evaluate the residual, given y and ydot.
61 /*!
62 * IDA allows passing in a void* pointer to access external data. Instead of
63 * requiring the user to provide a residual function directly to IDA (which
64 * would require using the sundials data types N_Vector, etc.), we define
65 * this function as the single function that IDA always calls. The real
66 * evaluation of the residual is done by an instance of a subclass of
67 * ResidEval, passed in to this function as a pointer in the parameters.
68 *
69 * FROM IDA WRITEUP -> What the IDA solver expects as a return flag from its
70 * residual routines:
71 *
72 * A IDAResFn res should return a value of 0 if successful, a positive value
73 * if a recoverable error occurred (for example, yy has an illegal value), or a
74 * negative value if a nonrecoverable error occurred. In the latter case, the
75 * program halts. If a recoverable error occurred, the integrator will
76 * attempt to correct and retry.
77 */
78 static int ida_resid(realtype t, N_Vector y, N_Vector ydot, N_Vector r, void* f_data)
79 {
81 Cantera::ResidJacEval* f = d->m_func;
82 Cantera::IDA_Solver* s = d->m_solver;
83 double delta_t = s->getCurrentStepFromIDA();
84 // TODO evaluate evalType. Assumed to be Base_ResidEval
85 int flag = f->evalResidNJ(t, delta_t, NV_DATA_S(y), NV_DATA_S(ydot),
86 NV_DATA_S(r));
87 if (flag < 0) {
88 // This signals to IDA that a nonrecoverable error has occurred.
89 return flag;
90 } else {
91 return 0;
92 }
93 }
94
95 //! Function called by by IDA to evaluate the Jacobian, given y and ydot.
96 /*!
97 * typedef int (*IDADlsDenseJacFn)(sd_size_t N, realtype t, realtype c_j,
98 * N_Vector y, N_Vector yp, N_Vector r,
99 * DlsMat Jac, void *user_data,
100 * N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
101 *
102 * A IDADlsDenseJacFn should return
103 * - 0 if successful,
104 * - a positive int if a recoverable error occurred, or
105 * - a negative int if a nonrecoverable error occurred.
106 *
107 * In the case of a recoverable error return, the integrator will attempt to
108 * recover by reducing the stepsize (which changes cj).
109 */
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)
114 {
116 Cantera::ResidJacEval* f = d->m_func;
117 Cantera::IDA_Solver* s = d->m_solver;
118 double delta_t = s->getCurrentStepFromIDA();
119 double** cols;
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);
124 } else {
125 return 1; // Unknown SUNMatrix type
126 }
127 f->evalJacobianDP(t, delta_t, c_j, NV_DATA_S(y), NV_DATA_S(yp),
128 cols, NV_DATA_S(r));
129 return 0;
130 }
131#else
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)
134 {
136 Cantera::ResidJacEval* f = d->m_func;
137 Cantera::IDA_Solver* s = d->m_solver;
138 double delta_t = s->getCurrentStepFromIDA();
139 f->evalJacobianDP(t, delta_t, c_j, NV_DATA_S(y), NV_DATA_S(ydot),
140 Jac->cols, NV_DATA_S(r));
141 return 0;
142 }
143#endif
144
145}
146
147namespace Cantera
148{
149
151 DAE_Solver(f),
152 m_ida_mem(0),
153 m_t0(0.0),
154 m_y(0),
155 m_ydot(0),
156 m_id(0),
157 m_constraints(0),
158 m_abstol(0),
159 m_type(0),
160 m_itol(IDA_SS),
161 m_iter(0),
162 m_reltol(1.e-9),
163 m_abstols(1.e-15),
164 m_nabs(0),
165 m_hmax(0.0),
166 m_hmin(0.0),
167 m_h0(0.0),
168 m_maxsteps(20000),
169 m_maxord(0),
170 m_formJac(0),
171 m_tstop(0.0),
172 m_told_old(0.0),
173 m_told(0.0),
174 m_tcurrent(0.0),
175 m_deltat(0.0),
176 m_maxErrTestFails(-1),
177 m_maxNonlinIters(0),
178 m_maxNonlinConvFails(-1),
179 m_setSuppressAlg(0),
180 m_mupper(0),
181 m_mlower(0)
182{
183}
184
185IDA_Solver::~IDA_Solver()
186{
187 if (m_ida_mem) {
188 IDAFree(&m_ida_mem);
189 }
190 if (m_y) {
191 N_VDestroy_Serial(m_y);
192 }
193 if (m_ydot) {
194 N_VDestroy_Serial(m_ydot);
195 }
196 if (m_abstol) {
197 N_VDestroy_Serial(m_abstol);
198 }
199 if (m_constraints) {
200 N_VDestroy_Serial(m_constraints);
201 }
202}
203
204doublereal IDA_Solver::solution(int k) const
205{
206 return NV_Ith_S(m_y,k);
207}
208
209const doublereal* IDA_Solver::solutionVector() const
210{
211 return NV_DATA_S(m_y);
212}
213
214doublereal IDA_Solver::derivative(int k) const
215{
216 return NV_Ith_S(m_ydot,k);
217}
218
219const doublereal* IDA_Solver::derivativeVector() const
220{
221 return NV_DATA_S(m_ydot);
222}
223
224void IDA_Solver::setTolerances(double reltol, double* abstol)
225{
226 m_itol = IDA_SV;
227 if (!m_abstol) {
228 #if CT_SUNDIALS_VERSION >= 60
229 m_abstol = N_VNew_Serial(m_neq, m_sundials_ctx.get());
230 #else
231 m_abstol = N_VNew_Serial(m_neq);
232 #endif
233 }
234 for (int i = 0; i < m_neq; i++) {
235 NV_Ith_S(m_abstol, i) = abstol[i];
236 }
237 m_reltol = reltol;
238 if (m_ida_mem) {
239 int flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
240 if (flag != IDA_SUCCESS) {
241 throw CanteraError("IDA_Solver::setTolerances",
242 "Memory allocation failed.");
243 }
244 }
245}
246
247void IDA_Solver::setTolerances(doublereal reltol, doublereal abstol)
248{
249 m_itol = IDA_SS;
250 m_reltol = reltol;
251 m_abstols = abstol;
252 if (m_ida_mem) {
253 int flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
254 if (flag != IDA_SUCCESS) {
255 throw CanteraError("IDA_Solver::setTolerances",
256 "Memory allocation failed.");
257 }
258 }
259}
260
261void IDA_Solver::setLinearSolverType(int solverType)
262{
263 m_type = solverType;
264}
265
267{
268 setLinearSolverType(0);
269}
270
271void IDA_Solver::setBandedLinearSolver(int m_upper, int m_lower)
272{
273 m_type = 2;
274 m_upper = m_mupper;
275 m_mlower = m_lower;
276}
277
278void IDA_Solver::setMaxOrder(int n)
279{
280 m_maxord = n;
281}
282
284{
285 m_maxsteps = n;
286}
287
289{
290 m_h0 = h0;
291}
292
293void IDA_Solver::setStopTime(doublereal tstop)
294{
295 m_tstop = tstop;
296}
297
299{
300 doublereal hcur;
301 IDAGetCurrentStep(m_ida_mem, &hcur);
302 return hcur;
303}
304
306{
307 m_formJac = formJac;
308 if (m_ida_mem && m_formJac == 1) {
309 #if CT_SUNDIALS_VERSION >= 60
310 int flag = IDASetJacFn(m_ida_mem, ida_jacobian);
311 #elif CT_SUNDIALS_VERSION >= 30
312 int flag = IDADlsSetJacFn(m_ida_mem, ida_jacobian);
313 #else
314 int flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
315 #endif
316 if (flag != IDA_SUCCESS) {
317 throw CanteraError("IDA_Solver::setJacobianType",
318 "IDADlsSetDenseJacFn failed.");
319 }
320 }
321}
322
323void IDA_Solver::setMaxErrTestFailures(int maxErrTestFails)
324{
325 m_maxErrTestFails = maxErrTestFails;
326}
327
329{
331}
332
334{
336}
337
338void IDA_Solver::inclAlgebraicInErrorTest(bool yesno)
339{
340 if (yesno) {
342 } else {
344 }
345}
346
347void IDA_Solver::init(doublereal t0)
348{
349 m_t0 = t0;
350 m_told = t0;
351 m_told_old = t0;
352 m_tcurrent = t0;
353 if (m_y) {
354 N_VDestroy_Serial(m_y);
355 }
356 if (m_ydot) {
357 N_VDestroy_Serial(m_ydot);
358 }
359 if (m_id) {
360 N_VDestroy_Serial(m_id);
361 }
362 if (m_constraints) {
363 N_VDestroy_Serial(m_constraints);
364 }
365
366 #if CT_SUNDIALS_VERSION >= 60
367 m_y = N_VNew_Serial(m_neq, m_sundials_ctx.get());
368 m_ydot = N_VNew_Serial(m_neq, m_sundials_ctx.get());
369 m_constraints = N_VNew_Serial(m_neq, m_sundials_ctx.get());
370 #else
371 m_y = N_VNew_Serial(m_neq);
372 m_ydot = N_VNew_Serial(m_neq);
373 m_constraints = N_VNew_Serial(m_neq);
374 #endif
375
376 for (int i=0; i<m_neq; i++) {
377 NV_Ith_S(m_y, i) = 0.0;
378 NV_Ith_S(m_ydot, i) = 0.0;
379 NV_Ith_S(m_constraints, i) = 0.0;
380 }
381
382 // get the initial conditions
383 m_resid.getInitialConditions(m_t0, NV_DATA_S(m_y), NV_DATA_S(m_ydot));
384
385 if (m_ida_mem) {
386 IDAFree(&m_ida_mem);
387 }
388
389 /* Call IDACreate */
390 #if CT_SUNDIALS_VERSION >= 60
391 m_ida_mem = IDACreate(m_sundials_ctx.get());
392 #else
393 m_ida_mem = IDACreate();
394 #endif
395
396 int flag = 0;
397 if (m_itol == IDA_SV) {
398 flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
399 if (flag != IDA_SUCCESS) {
400 if (flag == IDA_MEM_FAIL) {
401 throw CanteraError("IDA_Solver::init",
402 "Memory allocation failed.");
403 } else if (flag == IDA_ILL_INPUT) {
404 throw CanteraError("IDA_Solver::init",
405 "Illegal value for IDAInit input argument.");
406 } else {
407 throw CanteraError("IDA_Solver::init", "IDAInit failed.");
408 }
409 }
410 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
411 if (flag != IDA_SUCCESS) {
412 throw CanteraError("IDA_Solver::init", "Memory allocation failed.");
413 }
414 } else {
415 flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
416 if (flag != IDA_SUCCESS) {
417 if (flag == IDA_MEM_FAIL) {
418 throw CanteraError("IDA_Solver::init",
419 "Memory allocation failed.");
420 } else if (flag == IDA_ILL_INPUT) {
421 throw CanteraError("IDA_Solver::init",
422 "Illegal value for IDAInit input argument.");
423 } else {
424 throw CanteraError("IDA_Solver::init", "IDAInit failed.");
425 }
426 }
427 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
428 if (flag != IDA_SUCCESS) {
429 throw CanteraError("IDA_Solver::init", "Memory allocation failed.");
430 }
431 }
432
433 // set the linear solver type
434 if (m_type == 1 || m_type == 0) {
435 long int N = m_neq;
436 int flag;
437 #if CT_SUNDIALS_VERSION >= 30
438 SUNLinSolFree((SUNLinearSolver) m_linsol);
439 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
440 #if CT_SUNDIALS_VERSION >= 60
441 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
442 #else
443 m_linsol_matrix = SUNDenseMatrix(N, N);
444 #endif
445 if (m_linsol_matrix == nullptr) {
446 throw CanteraError("IDA_Solver::init",
447 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
448 }
449 #if CT_SUNDIALS_VERSION >= 60
450 #if CT_SUNDIALS_USE_LAPACK
451 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
452 m_sundials_ctx.get());
453 #else
454 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
455 m_sundials_ctx.get());
456 #endif
457 flag = IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
458 (SUNMatrix) m_linsol_matrix);
459 #else
460 #if CT_SUNDIALS_USE_LAPACK
461 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
462 #else
463 m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
464 #endif
465 flag = IDADlsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
466 (SUNMatrix) m_linsol_matrix);
467 #endif
468 #else
469 flag = IDADense(m_ida_mem, N);
470 #endif
471 if (flag) {
472 throw CanteraError("IDA_Solver::init", "IDADense failed");
473 }
474 } else if (m_type == 2) {
475 long int N = m_neq;
476 long int nu = m_mupper;
477 long int nl = m_mlower;
478 #if CT_SUNDIALS_VERSION >= 30
479 SUNLinSolFree((SUNLinearSolver) m_linsol);
480 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
481 #if CT_SUNDIALS_VERSION >= 60
482 m_linsol_matrix = SUNBandMatrix(N, nu, nl, m_sundials_ctx.get());
483 #elif CT_SUNDIALS_VERSION >= 40
484 m_linsol_matrix = SUNBandMatrix(N, nu, nl);
485 #else
486 m_linsol_matrix = SUNBandMatrix(N, nu, nl, nu+nl);
487 #endif
488 if (m_linsol_matrix == nullptr) {
489 throw CanteraError("IDA_Solver::init",
490 "Unable to create SUNBandMatrix of size {} with bandwidths "
491 "{} and {}", N, nu, nl);
492 }
493 #if CT_SUNDIALS_VERSION >= 60
494 #if CT_SUNDIALS_USE_LAPACK
495 m_linsol = SUNLinSol_LapackBand(m_y, (SUNMatrix) m_linsol_matrix,
496 m_sundials_ctx.get());
497 #else
498 m_linsol = SUNLinSol_Band(m_y, (SUNMatrix) m_linsol_matrix,
499 m_sundials_ctx.get());
500 #endif
501 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
502 (SUNMatrix) m_linsol_matrix);
503 #else
504 #if CT_SUNDIALS_USE_LAPACK
505 m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
506 #else
507 m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
508 #endif
509 IDADlsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
510 (SUNMatrix) m_linsol_matrix);
511 #endif
512 #else
513 IDABand(m_ida_mem, N, nu, nl);
514 #endif
515 } else {
516 throw CanteraError("IDA_Solver::init",
517 "unsupported linear solver type");
518 }
519
520 if (m_formJac == 1) {
521 #if CT_SUNDIALS_VERSION >= 60
522 flag = IDASetJacFn(m_ida_mem, ida_jacobian);
523 #elif CT_SUNDIALS_VERSION >= 30
524 flag = IDADlsSetJacFn(m_ida_mem, ida_jacobian);
525 #else
526 flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
527 #endif
528 if (flag != IDA_SUCCESS) {
529 throw CanteraError("IDA_Solver::init",
530 "IDADlsSetDenseJacFn failed.");
531 }
532 }
533
534 // pass a pointer to func in m_data
535 m_fdata.reset(new ResidData(&m_resid, this, m_resid.nparams()));
536 flag = IDASetUserData(m_ida_mem, m_fdata.get());
537 if (flag != IDA_SUCCESS) {
538 throw CanteraError("IDA_Solver::init", "IDASetUserData failed.");
539 }
540
541 // set options
542 if (m_maxord > 0) {
543 flag = IDASetMaxOrd(m_ida_mem, m_maxord);
544 if (flag != IDA_SUCCESS) {
545 throw CanteraError("IDA_Solver::init", "IDASetMaxOrd failed.");
546 }
547 }
548 if (m_maxsteps > 0) {
549 flag = IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
550 if (flag != IDA_SUCCESS) {
551 throw CanteraError("IDA_Solver::init", "IDASetMaxNumSteps failed.");
552 }
553 }
554 if (m_h0 > 0.0) {
555 flag = IDASetInitStep(m_ida_mem, m_h0);
556 if (flag != IDA_SUCCESS) {
557 throw CanteraError("IDA_Solver::init", "IDASetInitStep failed.");
558 }
559 }
560 if (m_tstop > 0.0) {
561 flag = IDASetStopTime(m_ida_mem, m_tstop);
562 if (flag != IDA_SUCCESS) {
563 throw CanteraError("IDA_Solver::init", "IDASetStopTime failed.");
564 }
565 }
566 if (m_maxErrTestFails >= 0) {
567 flag = IDASetMaxErrTestFails(m_ida_mem, m_maxErrTestFails);
568 if (flag != IDA_SUCCESS) {
569 throw CanteraError("IDA_Solver::init",
570 "IDASetMaxErrTestFails failed.");
571 }
572 }
573 if (m_maxNonlinIters > 0) {
574 flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
575 if (flag != IDA_SUCCESS) {
576 throw CanteraError("IDA_Solver::init",
577 "IDASetmaxNonlinIters failed.");
578 }
579 }
580 if (m_maxNonlinConvFails >= 0) {
581 flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
582 if (flag != IDA_SUCCESS) {
583 throw CanteraError("IDA_Solver::init",
584 "IDASetMaxConvFails failed.");
585 }
586 }
587 if (m_setSuppressAlg != 0) {
588 flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
589 if (flag != IDA_SUCCESS) {
590 throw CanteraError("IDA_Solver::init", "IDASetSuppressAlg failed.");
591 }
592 }
593}
594
595void IDA_Solver::correctInitial_Y_given_Yp(doublereal* y, doublereal* yp, doublereal tout)
596{
597 doublereal tout1 = tout;
598 if (tout == 0.0) {
599 double h0 = 1.0E-5;
600 if (m_h0 > 0.0) {
601 h0 = m_h0;
602 }
603 tout1 = m_t0 + h0;
604 }
605
606 int flag = IDACalcIC(m_ida_mem, IDA_Y_INIT, tout1);
607 if (flag != IDA_SUCCESS) {
608 throw CanteraError("IDA_Solver::correctInitial_Y_given_Yp",
609 "IDACalcIC failed: error = {}", flag);
610 }
611
612 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
613 if (flag != IDA_SUCCESS) {
614 throw CanteraError("IDA_Solver::correctInitial_Y_given_Yp",
615 "IDAGetSolution failed: error = {}", flag);
616 }
617 for (int i = 0; i < m_neq; i++) {
618 y[i] = NV_Ith_S(m_y, i);
619 yp[i] = NV_Ith_S(m_ydot, i);
620 }
621}
622
623void IDA_Solver::correctInitial_YaYp_given_Yd(doublereal* y, doublereal* yp, doublereal tout)
624{
625 int icopt = IDA_YA_YDP_INIT;
626 doublereal tout1 = tout;
627 if (tout == 0.0) {
628 double h0 = 1.0E-5;
629 if (m_h0 > 0.0) {
630 h0 = m_h0;
631 }
632 tout1 = m_t0 + h0;
633 }
634
635 int flag = IDACalcIC(m_ida_mem, icopt, tout1);
636 if (flag != IDA_SUCCESS) {
637 throw CanteraError("IDA_Solver::correctInitial_YaYp_given_Yd",
638 "IDACalcIC failed: error = {}", flag);
639 }
640
641 flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
642 if (flag != IDA_SUCCESS) {
643 throw CanteraError("IDA_Solver::correctInitial_YaYp_given_Yd",
644 "IDAGetSolution failed: error = {}", flag);
645 }
646 for (int i = 0; i < m_neq; i++) {
647 y[i] = NV_Ith_S(m_y, i);
648 yp[i] = NV_Ith_S(m_ydot, i);
649 }
650}
651
652int IDA_Solver::solve(double tout)
653{
654 double tretn = tout - 1000;
655 int flag;
656 flag = IDASetStopTime(m_ida_mem, tout);
657 if (flag != IDA_SUCCESS) {
658 throw CanteraError("IDA_Solver::solve", "IDA error encountered.");
659 }
660 while (tretn < tout) {
661 if (tout <= m_tcurrent) {
662 throw CanteraError("IDA_Solver::solve", "tout <= tcurrent");
663 }
666 flag = IDASolve(m_ida_mem, tout, &tretn, m_y, m_ydot, IDA_ONE_STEP);
667 if (flag < 0) {
668 throw CanteraError("IDA_Solver::solve", "IDA error encountered.");
669 } else if (flag == IDA_TSTOP_RETURN) {
670 // we've reached our goal, and have actually integrated past it
671 } else if (flag == IDA_ROOT_RETURN) {
672 // not sure what to do with this yet
673 } else if (flag == IDA_WARNING) {
674 throw CanteraError("IDA_Solver::solve", "IDA Warning encountered.");
675 }
676 m_tcurrent = tretn;
678 };
679
680 if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
681 throw CanteraError("IDA_Solver::solve", "IDA error encountered.");
682 }
683 return flag;
684}
685
686double IDA_Solver::step(double tout)
687{
688 double t;
689 int flag;
690 if (tout <= m_tcurrent) {
691 throw CanteraError("IDA_Solver::step", "tout <= tcurrent");
692 }
695 flag = IDASolve(m_ida_mem, tout, &t, m_y, m_ydot, IDA_ONE_STEP);
696 if (flag < 0) {
697 throw CanteraError("IDA_Solver::step", "IDA error encountered.");
698 } else if (flag == IDA_TSTOP_RETURN) {
699 // we've reached our goal, and have actually integrated past it
700 } else if (flag == IDA_ROOT_RETURN) {
701 // not sure what to do with this yet
702 } else if (flag == IDA_WARNING) {
703 throw CanteraError("IDA_Solver::step", "IDA Warning encountered.");
704 }
705 m_tcurrent = t;
707 return t;
708}
709
710doublereal IDA_Solver::getOutputParameter(int flag) const
711{
712 long int lenrw, leniw;
713 switch (flag) {
714 case REAL_WORKSPACE_SIZE:
715 flag = IDAGetWorkSpace(m_ida_mem, &lenrw, &leniw);
716 return doublereal(lenrw);
717 break;
718 }
719 return 0.0;
720}
721
722}
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.
Definition: IDA_Solver.cpp:132
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.
Definition: IDA_Solver.cpp:78
Header file for class IDA_Solver.
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
Wrapper for DAE solvers.
Definition: DAE_Solver.h:76
integer m_neq
Number of total equations in the system.
Definition: DAE_Solver.h:249
Wrapper for Sundials IDA solver.
Definition: IDA_Solver.h:42
virtual doublereal step(doublereal tout)
Take one internal step.
Definition: IDA_Solver.cpp:686
virtual void setStopTime(doublereal tstop)
Set the stop time.
Definition: IDA_Solver.cpp:293
virtual void correctInitial_Y_given_Yp(doublereal *y, doublereal *yp, doublereal tout)
Calculate consistent value of the starting solution given the starting solution derivatives.
Definition: IDA_Solver.cpp:595
virtual void setMaxNumSteps(int n)
Set the maximum number of time steps.
Definition: IDA_Solver.cpp:283
virtual void setBandedLinearSolver(int m_upper, int m_lower)
Set up the problem to use a band solver.
Definition: IDA_Solver.cpp:271
virtual doublereal solution(int k) const
the current value of solution component k.
Definition: IDA_Solver.cpp:204
N_Vector m_y
Current value of the solution vector.
Definition: IDA_Solver.h:235
virtual void setJacobianType(int formJac)
Set the form of the Jacobian.
Definition: IDA_Solver.cpp:305
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.
Definition: IDA_Solver.cpp:623
void * m_linsol
Sundials linear solver object.
Definition: IDA_Solver.h:227
int m_maxNonlinConvFails
Maximum number of nonlinear convergence failures.
Definition: IDA_Solver.h:299
virtual void setInitialStepSize(doublereal h0)
Set the initial step size.
Definition: IDA_Solver.cpp:288
virtual int solve(doublereal tout)
Step the system to a final value of the time.
Definition: IDA_Solver.cpp:652
int m_maxNonlinIters
Maximum number of nonlinear solver iterations at one solution.
Definition: IDA_Solver.h:296
SundialsContext m_sundials_ctx
SUNContext object for Sundials>=6.0.
Definition: IDA_Solver.h:229
IDA_Solver(ResidJacEval &f)
Constructor.
Definition: IDA_Solver.cpp:150
doublereal m_told
Value of the previous time.
Definition: IDA_Solver.h:281
virtual void setDenseLinearSolver()
Set up the problem to use a dense linear direct solver.
Definition: IDA_Solver.cpp:266
doublereal m_deltat
Value of deltaT for the current step.
Definition: IDA_Solver.h:287
int m_maxErrTestFails
maximum number of error test failures
Definition: IDA_Solver.h:290
virtual void setTolerances(doublereal reltol, doublereal *abstol)
Set error tolerances.
Definition: IDA_Solver.cpp:224
virtual void setMaxNonlinIterations(int n)
Set the maximum number of nonlinear iterations on a timestep.
Definition: IDA_Solver.cpp:328
doublereal m_tcurrent
Value of the current time.
Definition: IDA_Solver.h:284
int m_maxord
maximum time step order of the method
Definition: IDA_Solver.h:263
virtual doublereal getOutputParameter(int flag) const
Get the value of a solver-specific output parameter.
Definition: IDA_Solver.cpp:710
virtual void init(doublereal t0)
initialize.
Definition: IDA_Solver.cpp:347
int m_formJac
Form of the Jacobian.
Definition: IDA_Solver.h:272
doublereal m_told_old
Value of the previous, previous time.
Definition: IDA_Solver.h:278
virtual void setMaxNonlinConvFailures(int n)
Set the maximum number of nonlinear solver convergence failures.
Definition: IDA_Solver.cpp:333
virtual double getCurrentStepFromIDA()
Get the current step size from IDA via a call.
Definition: IDA_Solver.cpp:298
virtual doublereal derivative(int k) const
the current value of the derivative of solution component k.
Definition: IDA_Solver.cpp:214
doublereal m_h0
Value of the initial time step.
Definition: IDA_Solver.h:257
void * m_linsol_matrix
matrix used by Sundials
Definition: IDA_Solver.h:228
N_Vector m_ydot
Current value of the derivative of the solution vector.
Definition: IDA_Solver.h:238
doublereal m_t0
Initial value of the time.
Definition: IDA_Solver.h:232
int m_maxsteps
Maximum number of time steps allowed.
Definition: IDA_Solver.h:260
doublereal m_tstop
maximum time
Definition: IDA_Solver.h:275
int m_setSuppressAlg
If true, the algebraic variables don't contribute to error tolerances.
Definition: IDA_Solver.h:302
void * m_ida_mem
Pointer to the IDA memory for the problem.
Definition: IDA_Solver.h:226
A simple class to hold an array of parameter values and a pointer to an instance of a subclass of Res...
Definition: IDA_Solver.cpp:44
int nparams() const
Return the number of parameters in the calculation.
Definition: ResidEval.h:156
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
Definition: ResidJacEval.h:56
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.
Definition: AnyMap.h:29
Contains declarations for string manipulation functions within Cantera.