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 to return stats
181 AnyMap stats;
182
183 // long int linear solver stats provided by IDAS
184 long int steps = 0, stepSolveFails = 0, resEvals = 0, errTestFails = 0,
185 jacEvals = 0, linSetup = 0, linResEvals = 0, linIters = 0,
186 linConvFails = 0, precEvals = 0, precSolves = 0, jtSetupEvals = 0,
187 jTimesEvals = 0, nonlinIters = 0, nonlinConvFails = 0;
188 int lastOrder = 0;
189
190 int flag = IDAGetNumSteps(m_ida_mem, &steps);
191 checkError(flag, "solverStats", "IDAGetNumSteps");
192
193 // Remaining stats are best-effort; leave corresponding values at zero if the
194 // selected linear solver does not report a given counter.
195 IDAGetNumStepSolveFails(m_ida_mem, &stepSolveFails);
196 IDAGetNumResEvals(m_ida_mem, &resEvals);
197 IDAGetNumNonlinSolvIters(m_ida_mem, &nonlinIters);
198 IDAGetNumNonlinSolvConvFails(m_ida_mem, &nonlinConvFails);
199 IDAGetNumErrTestFails(m_ida_mem, &errTestFails);
200 IDAGetLastOrder(m_ida_mem, &lastOrder);
201
202 IDAGetNumJacEvals(m_ida_mem, &jacEvals);
203 IDAGetNumLinResEvals(m_ida_mem, &linResEvals);
204 IDAGetNumLinSolvSetups(m_ida_mem, &linSetup);
205 IDAGetNumLinIters(m_ida_mem, &linIters);
206 IDAGetNumLinConvFails(m_ida_mem, &linConvFails);
207 IDAGetNumPrecEvals(m_ida_mem, &precEvals);
208 IDAGetNumPrecSolves(m_ida_mem, &precSolves);
209 IDAGetNumJTSetupEvals(m_ida_mem, &jtSetupEvals);
210 IDAGetNumJtimesEvals(m_ida_mem, &jTimesEvals);
211
212 stats["steps"] = steps;
213 stats["step_solve_fails"] = stepSolveFails;
214 stats["res_evals"] = resEvals;
215 stats["rhs_evals"] = resEvals;
216 stats["nonlinear_iters"] = nonlinIters;
217 stats["nonlinear_conv_fails"] = nonlinConvFails;
218 stats["err_test_fails"] = errTestFails;
219 stats["last_order"] = lastOrder;
220
221 stats["jac_evals"] = jacEvals;
222 stats["lin_solve_setups"] = linSetup;
223 stats["lin_rhs_evals"] = linResEvals;
224 stats["lin_iters"] = linIters;
225 stats["lin_conv_fails"] = linConvFails;
226 stats["prec_evals"] = precEvals;
227 stats["prec_solves"] = precSolves;
228 stats["jt_vec_setup_evals"] = jtSetupEvals;
229 stats["jt_vec_prod_evals"] = jTimesEvals;
230 return stats;
231}
232
233void IdasIntegrator::setMaxNonlinIterations(int n)
234{
236 if (m_ida_mem) {
237 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
238 checkError(flag, "setMaxNonlinIterations", "IDASetMaxNonlinIters");
239 }
240}
241
242void IdasIntegrator::setMaxNonlinConvFailures(int n)
243{
245 if (m_ida_mem) {
246 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
247 checkError(flag, "setMaxNonlinConvFailures", "IDASetMaxConvFails");
248 }
249}
250
251void IdasIntegrator::includeAlgebraicInErrorTest(bool yesno)
252{
253 m_setSuppressAlg = !yesno;
254 if (m_ida_mem) {
255 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
256 checkError(flag, "inclAlgebraicInErrorTest", "IDASetSuppressAlg");
257 }
258}
259
261{
262 m_neq = func.neq();
263 m_t0 = t0;
264 m_time = t0;
265 m_tInteg = t0;
266 m_func = &func;
267 func.clearErrors();
268
269 if (m_y) {
270 N_VDestroy_Serial(m_y); // free solution vector if already allocated
271 }
272 m_y = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
273 N_VConst(0.0, m_y);
274
275 if (m_ydot) {
276 N_VDestroy_Serial(m_ydot); // free derivative vector if already allocated
277 }
278 m_ydot = newNVector(m_neq, m_sundials_ctx);
279 N_VConst(0.0, m_ydot);
280
281 // check abs tolerance array size
282 if (m_itol == IDA_SV && m_nabs < m_neq) {
283 throw CanteraError("IdasIntegrator::initialize",
284 "not enough absolute tolerance values specified.");
285 }
286
287 if (m_constraints) {
288 N_VDestroy_Serial(m_constraints);
289 }
290 m_constraints = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
291 // set the constraints
292 func.getConstraints(asSpan(m_constraints));
293
294 // get the initial conditions
296
297 if (m_ida_mem) {
298 IDAFree(&m_ida_mem);
299 }
300
301 //! Create the IDA solver
302 m_ida_mem = IDACreate(m_sundials_ctx.get());
303 if (!m_ida_mem) {
304 throw CanteraError("IdasIntegrator::initialize", "IDACreate failed.");
305 }
306
307 int flag = IDAInit(m_ida_mem, ida_rhs, m_t0, m_y, m_ydot);
308 if (flag != IDA_SUCCESS) {
309 if (flag == IDA_MEM_FAIL) {
310 throw CanteraError("IdasIntegrator::initialize",
311 "Memory allocation failed.");
312 } else if (flag == IDA_ILL_INPUT) {
313 throw CanteraError("IdasIntegrator::initialize",
314 "Illegal value for IDAInit input argument.");
315 } else {
316 throw CanteraError("IdasIntegrator::initialize", "IDAInit failed.");
317 }
318 }
319
320 #if SUNDIALS_VERSION_MAJOR >= 7
321 flag = SUNContext_PushErrHandler(m_sundials_ctx.get(), &sundials_err, this);
322 #else
323 flag = IDASetErrHandlerFn(m_ida_mem, &ida_err, this);
324 #endif
325
326 // set constraints
327 flag = IDASetId(m_ida_mem, m_constraints);
328 checkError(flag, "initialize", "IDASetId");
329
330 if (m_itol == IDA_SV) {
331 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
332 checkError(flag, "initialize", "IDASVtolerances");
333 } else {
334 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
335 checkError(flag, "initialize", "IDASStolerances");
336 }
337
338 flag = IDASetUserData(m_ida_mem, &func);
339 checkError(flag, "initialize", "IDASetUserData");
340
341 if (func.nparams() > 0) {
342 throw CanteraError("IdasIntegrator::initialize", "Sensitivity analysis "
343 "for DAE systems is not fully implemented");
344 sensInit(t0, func);
345 flag = IDASetSensParams(m_ida_mem, func.m_sens_params.data(),
346 func.m_paramScales.data(), NULL);
347 checkError(flag, "initialize", "IDASetSensParams");
348 }
349 applyOptions();
350}
351
352void IdasIntegrator::reinitialize(double t0, FuncEval& func)
353{
354 m_t0 = t0;
355 m_time = t0;
356 m_tInteg = t0;
358 m_func = &func;
359 func.clearErrors();
360
361 int result = IDAReInit(m_ida_mem, m_t0, m_y, m_ydot);
362 checkError(result, "reinitialize", "IDAReInit");
363 applyOptions();
364}
365
367{
368 if (m_type == "DENSE") {
369 sd_size_t N = static_cast<sd_size_t>(m_neq);
370 SUNLinSolFree((SUNLinearSolver) m_linsol);
371 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
372 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
373 #if CT_SUNDIALS_USE_LAPACK
374 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
375 m_sundials_ctx.get());
376 #else
377 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
378 m_sundials_ctx.get());
379 #endif
380 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
381 (SUNMatrix) m_linsol_matrix);
382 } else if (m_type == "GMRES") {
383 SUNLinSolFree((SUNLinearSolver) m_linsol);
384 m_linsol = SUNLinSol_SPGMR(m_y, SUN_PREC_NONE, 0, m_sundials_ctx.get());
385 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol, nullptr);
386 } else {
387 throw CanteraError("IdasIntegrator::applyOptions",
388 "unsupported linear solver flag '{}'", m_type);
389 }
390
391 if (m_init_step > 0) {
392 IDASetInitStep(m_ida_mem, m_init_step);
393 }
394
395 if (m_maxord > 0) {
396 int flag = IDASetMaxOrd(m_ida_mem, m_maxord);
397 checkError(flag, "applyOptions", "IDASetMaxOrd");
398 }
399 if (m_maxsteps > 0) {
400 IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
401 }
402 if (m_hmax > 0) {
403 IDASetMaxStep(m_ida_mem, m_hmax);
404 }
405 if (m_maxNonlinIters > 0) {
406 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
407 checkError(flag, "applyOptions", "IDASetMaxNonlinIters");
408 }
409 if (m_maxNonlinConvFails > 0) {
410 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
411 checkError(flag, "applyOptions", "IDASetMaxConvFails");
412 }
413 if (m_setSuppressAlg != 0) {
414 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
415 checkError(flag, "applyOptions", "IDASetSuppressAlg");
416 }
417}
418
419void IdasIntegrator::sensInit(double t0, FuncEval& func)
420{
421 m_np = func.nparams();
422 m_sens_ok = false;
423
424 N_Vector y = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
425 m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
426 for (size_t n = 0; n < m_np; n++) {
427 N_VConst(0.0, m_yS[n]);
428 }
429 N_VDestroy_Serial(y);
430 N_Vector ydot = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
431 m_ySdot = N_VCloneVectorArray(static_cast<int>(m_np), ydot);
432 for (size_t n = 0; n < m_np; n++) {
433 N_VConst(0.0, m_ySdot[n]);
434 }
435
436 int flag = IDASensInit(m_ida_mem, static_cast<sd_size_t>(m_np),
437 IDA_STAGGERED, IDASensResFn(0), m_yS, m_ySdot);
438 checkError(flag, "sensInit", "IDASensInit");
439
440 vector<double> atol(m_np);
441 for (size_t n = 0; n < m_np; n++) {
442 // This scaling factor is tuned so that reaction and species enthalpy
443 // sensitivities can be computed simultaneously with the same abstol.
444 atol[n] = m_abstolsens / func.m_paramScales[n];
445 }
446 flag = IDASensSStolerances(m_ida_mem, m_reltolsens, atol.data());
447 checkError(flag, "sensInit", "IDASensSStolerances");
448}
449
451{
452 if (tout == m_time) {
453 return;
454 } else if (tout < m_time) {
455 throw CanteraError("IdasIntegrator::integrate",
456 "Cannot integrate backwards in time.\n"
457 "Requested time {} < current time {}",
458 tout, m_time);
459 }
460 int nsteps = 0;
461 while (m_tInteg < tout) {
462 if (nsteps >= m_maxsteps) {
463 throw CanteraError("IdasIntegrator::integrate",
464 "Maximum number of timesteps ({}) taken without reaching output "
465 "time ({}).\nCurrent integrator time: {}",
466 nsteps, tout, m_time);
467 }
468 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
469 if (flag != IDA_SUCCESS) {
470 string f_errs = m_func->getErrors();
471 if (!f_errs.empty()) {
472 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
473 }
474 throw CanteraError("IdasIntegrator::integrate",
475 "IDA error encountered. Error code: {}\n{}\n"
476 "{}"
477 "Components with largest weighted error estimates:\n{}",
478 flag, m_error_message, f_errs, getErrorInfo(10));
479 }
480 nsteps++;
481 }
482 int flag = IDAGetDky(m_ida_mem, tout, 0, m_y);
483 checkError(flag, "integrate", "IDAGetDky");
484 m_time = tout;
485}
486
487double IdasIntegrator::step(double tout)
488{
489 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
490 if (flag != IDA_SUCCESS) {
491 string f_errs = m_func->getErrors();
492 if (!f_errs.empty()) {
493 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
494 }
495 throw CanteraError("IdasIntegrator::step",
496 "IDA error encountered. Error code: {}\n{}\n"
497 "{}"
498 "Components with largest weighted error estimates:\n{}",
499 flag, f_errs, m_error_message, getErrorInfo(10));
500
501 }
503 return m_time;
504}
505
506double IdasIntegrator::sensitivity(size_t k, size_t p)
507{
508 if (m_time == m_t0) {
509 // calls to IDAGetSens are only allowed after a successful time step.
510 return 0.0;
511 }
512 if (!m_sens_ok && m_np) {
513 int flag = IDAGetSensDky(m_ida_mem, m_time, 0, m_yS);
514 checkError(flag, "sensitivity", "IDAGetSens");
515 m_sens_ok = true;
516 }
517
518 if (k >= m_neq) {
519 throw CanteraError("IdasIntegrator::sensitivity",
520 "sensitivity: k out of range ({})", k);
521 }
522 if (p >= m_np) {
523 throw CanteraError("IdasIntegrator::sensitivity",
524 "sensitivity: p out of range ({})", p);
525 }
526 return NV_Ith_S(m_yS[p],k);
527}
528
530{
531 N_Vector errs = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
532 N_Vector errw = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
533 IDAGetErrWeights(m_ida_mem, errw);
534 IDAGetEstLocalErrors(m_ida_mem, errs);
535
536 vector<tuple<double, double, size_t>> weightedErrors;
537 for (size_t i=0; i<m_neq; i++) {
538 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
539 weightedErrors.emplace_back(-abs(err), err, i);
540 }
541 N_VDestroy(errs);
542 N_VDestroy(errw);
543
544 N = std::min(N, static_cast<int>(m_neq));
545 sort(weightedErrors.begin(), weightedErrors.end());
546 fmt::memory_buffer s;
547 for (int i=0; i<N; i++) {
548 fmt_append(s, "{}: {}\n", get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
549 }
550 return to_string(s);
551}
552
553void IdasIntegrator::checkError(long flag, const string& ctMethod,
554 const string& idaMethod) const
555{
556 if (flag == IDA_SUCCESS) {
557 return;
558 } else if (flag == IDA_MEM_NULL) {
559 throw CanteraError("IdasIntegrator::" + ctMethod,
560 "IDAS integrator is not initialized");
561 } else {
562 const char* flagname = IDAGetReturnFlagName(flag);
563 throw CanteraError("IdasIntegrator::" + ctMethod,
564 "{} returned error code {} ({}):\n{}",
565 idaMethod, flag, flagname, m_error_message);
566 }
567}
568
570{
571 if (t != BDF_Method) {
572 // IDA only has the BDF method
573 throw CanteraError("IdasIntegrator::setMethod", "unknown method");
574 }
575}
576
577}
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.