IdasIntegrator.cpp Source File#

Cantera: IdasIntegrator.cpp Source File
IdasIntegrator.cpp
Go to the documentation of this file.
1//! @file IdasIntegrator.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at http://www.cantera.org/license.txt for license and copyright information.
5
8
9#include "cantera/numerics/sundials_headers.h"
10
11using namespace std;
12
13namespace {
14
15N_Vector newNVector(size_t N, Cantera::SundialsContext& context)
16{
17#if CT_SUNDIALS_VERSION >= 60
18 return N_VNew_Serial(static_cast<sd_size_t>(N), context.get());
19#else
20 return N_VNew_Serial(static_cast<sd_size_t>(N));
21#endif
22}
23
24} // end anonymous namespace
25
26namespace Cantera
27{
28
29extern "C" {
30
31//! Function called by IDA to evaluate the residual, given y and ydot.
32/*!
33* IDA allows passing in a void* pointer to access external data. Instead of requiring
34* the user to provide a residual function directly to IDA (which would require using the
35* Sundials data types N_Vector, etc.), we define this function as the single function
36* that IDA always calls. The real evaluation of the residual is done by the
37* FuncEval::evalDae() method of an instance of a subclass of FuncEval that is passed
38* into this function as the `f_data` parameter.
39*/
40static int ida_rhs(realtype t, N_Vector y, N_Vector ydot, N_Vector r, void* f_data)
41{
42 FuncEval* f = (FuncEval*) f_data;
43 return f->evalDaeNoThrow(t, NV_DATA_S(y), NV_DATA_S(ydot), NV_DATA_S(r));
44}
45
46//! Function called by IDA when an error is encountered instead of writing to stdout.
47//! Here, save the error message provided by IDA so that it can be included in the
48//! subsequently raised CanteraError.
49static void ida_err(int error_code, const char* module,
50 const char* function, char* msg, void* eh_data)
51{
52 IdasIntegrator* integrator = (IdasIntegrator*) eh_data;
53 integrator->m_error_message = msg;
54 integrator->m_error_message += "\n";
55}
56
57}
58
60 : m_itol(IDA_SS)
61{
62}
63
64IdasIntegrator::~IdasIntegrator()
65{
66 if (m_ida_mem) {
67 IDAFree(&m_ida_mem);
68 }
69 if (m_y) {
70 N_VDestroy_Serial(m_y);
71 }
72 if (m_ydot) {
73 N_VDestroy_Serial(m_ydot);
74 }
75 if (m_abstol) {
76 N_VDestroy_Serial(m_abstol);
77 }
78 if (m_constraints) {
79 N_VDestroy_Serial(m_constraints);
80 }
81}
82
83double& IdasIntegrator::solution(size_t k)
84{
85 return NV_Ith_S(m_y, k);
86}
87
89{
90 return NV_DATA_S(m_y);
91}
92
93void IdasIntegrator::setTolerances(double reltol, size_t n, double* abstol)
94{
95 m_itol = IDA_SV;
96 if (n != m_nabs) {
97 m_nabs = n;
98 if (m_abstol) {
99 N_VDestroy_Serial(m_abstol);
100 }
101 m_abstol = newNVector(static_cast<sd_size_t>(n), m_sundials_ctx);
102 }
103 for (size_t i=0; i<n; i++) {
104 NV_Ith_S(m_abstol, i) = abstol[i];
105 }
106 m_reltol = reltol;
107}
108
109void IdasIntegrator::setTolerances(double reltol, double abstol)
110{
111 m_itol = IDA_SS;
112 m_reltol = reltol;
113 m_abstols = abstol;
114}
115
116void IdasIntegrator::setSensitivityTolerances(double reltol, double abstol)
117{
118 m_reltolsens = reltol;
119 m_abstolsens = abstol;
120}
121
122
123void IdasIntegrator::setLinearSolverType(const string& linearSolverType)
124{
126}
127
129{
130 if (m_ida_mem) {
131 int flag = IDASetMaxOrd(m_ida_mem, n);
132 checkError(flag, "setMaxOrder", "IDASetMaxOrd");
133 }
134 m_maxord = n;
135}
136
138{
139 m_hmax = hmax;
140 if (m_ida_mem) {
141 int flag = IDASetMaxStep(m_ida_mem, hmax);
142 checkError(flag, "setMaxStepSize", "IDASetMaxStep");
143 }
144}
145
147{
148 m_maxsteps = nmax;
149 if (m_ida_mem) {
150 IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
151 }
152}
153
155{
156 return m_maxsteps;
157}
158
160{
162 if (m_ida_mem) {
163 IDASetMaxErrTestFails(m_ida_mem, n);
164 }
165}
166
168{
169 AnyMap stats;
170 long int val;
171 int lastOrder;
172
173 int flag = IDAGetNumSteps(m_ida_mem, &val);
174 checkError(flag, "solverStats", "IDAGetNumSteps");
175 stats["steps"] = val;
176 IDAGetNumResEvals(m_ida_mem, &val);
177 stats["res_evals"] = val;
178 IDAGetNumLinSolvSetups(m_ida_mem, &val);
179 stats["lin_solve_setups"] = val;
180 IDAGetNumErrTestFails(m_ida_mem, &val);
181 stats["err_tests_fails"] = val;
182 IDAGetLastOrder(m_ida_mem, &lastOrder);
183 stats["last_order"] = lastOrder;
184 IDAGetNumNonlinSolvIters(m_ida_mem, &val);
185 stats["nonlinear_iters"] = val;
186 IDAGetNumNonlinSolvConvFails(m_ida_mem, &val);
187 stats["nonlinear_conv_fails"] = val;
188 return stats;
189}
190
191void IdasIntegrator::setMaxNonlinIterations(int n)
192{
194 if (m_ida_mem) {
195 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
196 checkError(flag, "setMaxNonlinIterations", "IDASetMaxNonlinIters");
197 }
198}
199
200void IdasIntegrator::setMaxNonlinConvFailures(int n)
201{
203 if (m_ida_mem) {
204 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
205 checkError(flag, "setMaxNonlinConvFailures", "IDASetMaxConvFails");
206 }
207}
208
209void IdasIntegrator::includeAlgebraicInErrorTest(bool yesno)
210{
211 m_setSuppressAlg = !yesno;
212 if (m_ida_mem) {
213 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
214 checkError(flag, "inclAlgebraicInErrorTest", "IDASetSuppressAlg");
215 }
216}
217
219{
220 m_neq = func.neq();
221 m_t0 = t0;
222 m_time = t0;
223 m_tInteg = t0;
224 m_func = &func;
225 func.clearErrors();
226
227 if (m_y) {
228 N_VDestroy_Serial(m_y); // free solution vector if already allocated
229 }
230 m_y = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
231 N_VConst(0.0, m_y);
232
233 if (m_ydot) {
234 N_VDestroy_Serial(m_ydot); // free derivative vector if already allocated
235 }
236 m_ydot = newNVector(m_neq, m_sundials_ctx);
237 N_VConst(0.0, m_ydot);
238
239 // check abs tolerance array size
240 if (m_itol == IDA_SV && m_nabs < m_neq) {
241 throw CanteraError("IdasIntegrator::initialize",
242 "not enough absolute tolerance values specified.");
243 }
244
245 if (m_constraints) {
246 N_VDestroy_Serial(m_constraints);
247 }
248 m_constraints = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
249 // set the constraints
250 func.getConstraints(NV_DATA_S(m_constraints));
251
252 // get the initial conditions
253 func.getStateDae(NV_DATA_S(m_y), NV_DATA_S(m_ydot));
254
255 if (m_ida_mem) {
256 IDAFree(&m_ida_mem);
257 }
258
259 //! Create the IDA solver
260 #if CT_SUNDIALS_VERSION >= 60
261 m_ida_mem = IDACreate(m_sundials_ctx.get());
262 #else
263 m_ida_mem = IDACreate();
264 #endif
265 if (!m_ida_mem) {
266 throw CanteraError("IdasIntegrator::initialize", "IDACreate failed.");
267 }
268
269 int flag = IDAInit(m_ida_mem, ida_rhs, m_t0, m_y, m_ydot);
270 if (flag != IDA_SUCCESS) {
271 if (flag == IDA_MEM_FAIL) {
272 throw CanteraError("IdasIntegrator::initialize",
273 "Memory allocation failed.");
274 } else if (flag == IDA_ILL_INPUT) {
275 throw CanteraError("IdasIntegrator::initialize",
276 "Illegal value for IDAInit input argument.");
277 } else {
278 throw CanteraError("IdasIntegrator::initialize", "IDAInit failed.");
279 }
280 }
281
282 flag = IDASetErrHandlerFn(m_ida_mem, &ida_err, this);
283
284 // set constraints
285 flag = IDASetId(m_ida_mem, m_constraints);
286 checkError(flag, "initialize", "IDASetId");
287
288 if (m_itol == IDA_SV) {
289 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
290 checkError(flag, "initialize", "IDASVtolerances");
291 } else {
292 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
293 checkError(flag, "initialize", "IDASStolerances");
294 }
295
296 flag = IDASetUserData(m_ida_mem, &func);
297 checkError(flag, "initialize", "IDASetUserData");
298
299 if (func.nparams() > 0) {
300 throw CanteraError("IdasIntegrator::initialize", "Sensitivity analysis "
301 "for DAE systems is not fully implemented");
302 sensInit(t0, func);
303 flag = IDASetSensParams(m_ida_mem, func.m_sens_params.data(),
304 func.m_paramScales.data(), NULL);
305 checkError(flag, "initialize", "IDASetSensParams");
306 }
307 applyOptions();
308}
309
310void IdasIntegrator::reinitialize(double t0, FuncEval& func)
311{
312 m_t0 = t0;
313 m_time = t0;
314 m_tInteg = t0;
315 func.getStateDae(NV_DATA_S(m_y), NV_DATA_S(m_ydot));
316 m_func = &func;
317 func.clearErrors();
318
319 int result = IDAReInit(m_ida_mem, m_t0, m_y, m_ydot);
320 checkError(result, "reinitialize", "IDAReInit");
321 applyOptions();
322}
323
325{
326 if (m_type == "DENSE") {
327 sd_size_t N = static_cast<sd_size_t>(m_neq);
328 SUNLinSolFree((SUNLinearSolver) m_linsol);
329 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
330 #if CT_SUNDIALS_VERSION >= 60
331 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
332 #else
333 m_linsol_matrix = SUNDenseMatrix(N, N);
334 #endif
335 #if CT_SUNDIALS_VERSION >= 60
336 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
337 #else
338 m_linsol_matrix = SUNDenseMatrix(N, N);
339 #endif
340 #if CT_SUNDIALS_USE_LAPACK
341 #if CT_SUNDIALS_VERSION >= 60
342 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
343 m_sundials_ctx.get());
344 #else
345 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
346 #endif
347 #else
348 #if CT_SUNDIALS_VERSION >= 60
349 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
350 m_sundials_ctx.get());
351 #else
352 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix);
353 #endif
354 #endif
355 #if CT_SUNDIALS_VERSION >= 60
356 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
357 (SUNMatrix) m_linsol_matrix);
358 #else
359 IDADlsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
360 (SUNMatrix) m_linsol_matrix);
361 #endif
362 } else if (m_type == "GMRES") {
363 #if CT_SUNDIALS_VERSION >= 60
364 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0, m_sundials_ctx.get());
365 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol, nullptr);
366 #elif CT_SUNDIALS_VERSION >= 40
367 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
368 IDASpilsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol);
369 #else
370 m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
371 IDASpilsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol);
372 #endif
373 } else {
374 throw CanteraError("IdasIntegrator::applyOptions",
375 "unsupported linear solver flag '{}'", m_type);
376 }
377
378 if (m_init_step > 0) {
379 IDASetInitStep(m_ida_mem, m_init_step);
380 }
381
382 if (m_maxord > 0) {
383 int flag = IDASetMaxOrd(m_ida_mem, m_maxord);
384 checkError(flag, "applyOptions", "IDASetMaxOrd");
385 }
386 if (m_maxsteps > 0) {
387 IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
388 }
389 if (m_hmax > 0) {
390 IDASetMaxStep(m_ida_mem, m_hmax);
391 }
392 if (m_maxNonlinIters > 0) {
393 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
394 checkError(flag, "applyOptions", "IDASetMaxNonlinIters");
395 }
396 if (m_maxNonlinConvFails > 0) {
397 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
398 checkError(flag, "applyOptions", "IDASetMaxConvFails");
399 }
400 if (m_setSuppressAlg != 0) {
401 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
402 checkError(flag, "applyOptions", "IDASetSuppressAlg");
403 }
404}
405
406void IdasIntegrator::sensInit(double t0, FuncEval& func)
407{
408 m_np = func.nparams();
409 m_sens_ok = false;
410
411 N_Vector y = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
412 #if CT_SUNDIALS_VERSION >= 60
413 m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
414 #else
415 m_yS = N_VCloneVectorArray_Serial(static_cast<int>(m_np), y);
416 #endif
417 for (size_t n = 0; n < m_np; n++) {
418 N_VConst(0.0, m_yS[n]);
419 }
420 N_VDestroy_Serial(y);
421 N_Vector ydot = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
422 #if CT_SUNDIALS_VERSION >= 60
423 m_ySdot = N_VCloneVectorArray(static_cast<int>(m_np), ydot);
424 #else
425 m_ySdot = N_VCloneVectorArray_Serial(static_cast<int>(m_np), ydot);
426 #endif
427 for (size_t n = 0; n < m_np; n++) {
428 N_VConst(0.0, m_ySdot[n]);
429 }
430
431 int flag = IDASensInit(m_ida_mem, static_cast<sd_size_t>(m_np),
432 IDA_STAGGERED, IDASensResFn(0), m_yS, m_ySdot);
433 checkError(flag, "sensInit", "IDASensInit");
434
435 vector<double> atol(m_np);
436 for (size_t n = 0; n < m_np; n++) {
437 // This scaling factor is tuned so that reaction and species enthalpy
438 // sensitivities can be computed simultaneously with the same abstol.
439 atol[n] = m_abstolsens / func.m_paramScales[n];
440 }
441 flag = IDASensSStolerances(m_ida_mem, m_reltolsens, atol.data());
442 checkError(flag, "sensInit", "IDASensSStolerances");
443}
444
446{
447 if (tout == m_time) {
448 return;
449 } else if (tout < m_time) {
450 throw CanteraError("IdasIntegrator::integrate",
451 "Cannot integrate backwards in time.\n"
452 "Requested time {} < current time {}",
453 tout, m_time);
454 }
455 int nsteps = 0;
456 while (m_tInteg < tout) {
457 if (nsteps >= m_maxsteps) {
458 throw CanteraError("IdasIntegrator::integrate",
459 "Maximum number of timesteps ({}) taken without reaching output "
460 "time ({}).\nCurrent integrator time: {}",
461 nsteps, tout, m_time);
462 }
463 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
464 if (flag != IDA_SUCCESS) {
465 string f_errs = m_func->getErrors();
466 if (!f_errs.empty()) {
467 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
468 }
469 throw CanteraError("IdasIntegrator::integrate",
470 "IDA error encountered. Error code: {}\n{}\n"
471 "{}"
472 "Components with largest weighted error estimates:\n{}",
473 flag, m_error_message, f_errs, getErrorInfo(10));
474 }
475 nsteps++;
476 }
477 int flag = IDAGetDky(m_ida_mem, tout, 0, m_y);
478 checkError(flag, "integrate", "IDAGetDky");
479 m_time = tout;
480}
481
482double IdasIntegrator::step(double tout)
483{
484 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
485 if (flag != IDA_SUCCESS) {
486 string f_errs = m_func->getErrors();
487 if (!f_errs.empty()) {
488 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
489 }
490 throw CanteraError("IdasIntegrator::step",
491 "IDA error encountered. Error code: {}\n{}\n"
492 "{}"
493 "Components with largest weighted error estimates:\n{}",
494 flag, f_errs, m_error_message, getErrorInfo(10));
495
496 }
498 return m_time;
499}
500
501double IdasIntegrator::sensitivity(size_t k, size_t p)
502{
503 if (m_time == m_t0) {
504 // calls to IDAGetSens are only allowed after a successful time step.
505 return 0.0;
506 }
507 if (!m_sens_ok && m_np) {
508 int flag = IDAGetSensDky(m_ida_mem, m_time, 0, m_yS);
509 checkError(flag, "sensitivity", "IDAGetSens");
510 m_sens_ok = true;
511 }
512
513 if (k >= m_neq) {
514 throw CanteraError("IdasIntegrator::sensitivity",
515 "sensitivity: k out of range ({})", k);
516 }
517 if (p >= m_np) {
518 throw CanteraError("IdasIntegrator::sensitivity",
519 "sensitivity: p out of range ({})", p);
520 }
521 return NV_Ith_S(m_yS[p],k);
522}
523
525{
526 N_Vector errs = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
527 N_Vector errw = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
528 IDAGetErrWeights(m_ida_mem, errw);
529 IDAGetEstLocalErrors(m_ida_mem, errs);
530
531 vector<tuple<double, double, size_t>> weightedErrors;
532 for (size_t i=0; i<m_neq; i++) {
533 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
534 weightedErrors.emplace_back(-abs(err), err, i);
535 }
536 N_VDestroy(errs);
537 N_VDestroy(errw);
538
539 N = std::min(N, static_cast<int>(m_neq));
540 sort(weightedErrors.begin(), weightedErrors.end());
541 fmt::memory_buffer s;
542 for (int i=0; i<N; i++) {
543 fmt_append(s, "{}: {}\n", get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
544 }
545 return to_string(s);
546}
547
548void IdasIntegrator::checkError(long flag, const string& ctMethod,
549 const string& idaMethod) const
550{
551 if (flag == IDA_SUCCESS) {
552 return;
553 } else if (flag == IDA_MEM_NULL) {
554 throw CanteraError("IdasIntegrator::" + ctMethod,
555 "IDAS integrator is not initialized");
556 } else {
557 const char* flagname = IDAGetReturnFlagName(flag);
558 throw CanteraError("IdasIntegrator::" + ctMethod,
559 "{} returned error code {} ({}):\n{}",
560 idaMethod, flag, flagname, m_error_message);
561 }
562}
563
565{
566 if (t != BDF_Method) {
567 // IDA only has the BDF method
568 throw CanteraError("IdasIntegrator::setMethod", "unknown method");
569 }
570}
571
572}
Header file for class IdasIntegrator.
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
Base class for exceptions thrown by Cantera classes.
Virtual base class for ODE/DAE right-hand-side function evaluators.
Definition FuncEval.h:32
virtual void getStateDae(double *y, double *ydot)
Fill in the vectors y and ydot with the current state of the system.
Definition FuncEval.h:148
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
Definition FuncEval.h:184
int evalDaeNoThrow(double t, double *y, double *ydot, double *residual)
Evaluate the right-hand side using return code to indicate status.
Definition FuncEval.cpp:39
void clearErrors()
Clear any previously-stored suppressed errors.
Definition FuncEval.h:174
virtual size_t neq() const =0
Number of equations.
virtual void getConstraints(double *constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Definition FuncEval.h:63
virtual size_t nparams() const
Number of sensitivity parameters.
Definition FuncEval.h:156
string getErrors() const
Return a string containing the text of any suppressed errors.
Definition FuncEval.cpp:71
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition FuncEval.h:181
Wrapper for Sundials IDAS solver.
double step(double tout) override
Integrate the system of equations.
void setMaxStepSize(double hmax) override
Set the maximum step size.
double * solution() override
The current value of the solution of the system of equations.
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. Zero means infinity.
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 setTolerances(double reltol, size_t n, double *abstol) override
Set error tolerances.
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 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.
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.
Definition Integrator.h:131
virtual int lastOrder() const
Order used during the last solution step.
Definition Integrator.h:188
A wrapper for managing a SUNContext object, need for Sundials >= 6.0.
void fmt_append(fmt::memory_buffer &b, Args... args)
Versions 6.2.0 and 6.2.1 of fmtlib do not include this define before they include windows....
Definition fmt.h:29
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
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.
Definition Integrator.h:23
@ BDF_Method
Backward Differentiation.
Definition Integrator.h:24
static int ida_rhs(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.
Contains declarations for string manipulation functions within Cantera.