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