Cantera  3.1.0b1
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#if SUNDIALS_VERSION_MAJOR >= 6
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(sunrealtype 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//! Function called by CVodes when an error is encountered instead of
58//! writing to stdout. Here, save the error message provided by CVodes so
59//! that it can be included in the subsequently raised CanteraError. Used by
60//! SUNDIALS 7.0 and newer.
61#if SUNDIALS_VERSION_MAJOR >= 7
62 static void sundials_err(int line, const char *func, const char *file,
63 const char *msg, SUNErrCode err_code,
64 void *err_user_data, SUNContext sunctx)
65 {
66 IdasIntegrator* integrator = (IdasIntegrator*) err_user_data;
67 integrator->m_error_message = fmt::format("{}: {}\n", func, msg);
68 }
69#endif
70
71}
72
74 : m_itol(IDA_SS)
75{
76}
77
78IdasIntegrator::~IdasIntegrator()
79{
80 if (m_ida_mem) {
81 IDAFree(&m_ida_mem);
82 }
83 if (m_y) {
84 N_VDestroy_Serial(m_y);
85 }
86 if (m_ydot) {
87 N_VDestroy_Serial(m_ydot);
88 }
89 if (m_abstol) {
90 N_VDestroy_Serial(m_abstol);
91 }
92 if (m_constraints) {
93 N_VDestroy_Serial(m_constraints);
94 }
95}
96
97double& IdasIntegrator::solution(size_t k)
98{
99 return NV_Ith_S(m_y, k);
100}
101
103{
104 return NV_DATA_S(m_y);
105}
106
107void IdasIntegrator::setTolerances(double reltol, size_t n, double* abstol)
108{
109 m_itol = IDA_SV;
110 if (n != m_nabs) {
111 m_nabs = n;
112 if (m_abstol) {
113 N_VDestroy_Serial(m_abstol);
114 }
115 m_abstol = newNVector(static_cast<sd_size_t>(n), m_sundials_ctx);
116 }
117 for (size_t i=0; i<n; i++) {
118 NV_Ith_S(m_abstol, i) = abstol[i];
119 }
120 m_reltol = reltol;
121}
122
123void IdasIntegrator::setTolerances(double reltol, double abstol)
124{
125 m_itol = IDA_SS;
126 m_reltol = reltol;
127 m_abstols = abstol;
128}
129
130void IdasIntegrator::setSensitivityTolerances(double reltol, double abstol)
131{
132 m_reltolsens = reltol;
133 m_abstolsens = abstol;
134}
135
136
137void IdasIntegrator::setLinearSolverType(const string& linearSolverType)
138{
140}
141
143{
144 if (m_ida_mem) {
145 int flag = IDASetMaxOrd(m_ida_mem, n);
146 checkError(flag, "setMaxOrder", "IDASetMaxOrd");
147 }
148 m_maxord = n;
149}
150
152{
153 m_hmax = hmax;
154 if (m_ida_mem) {
155 int flag = IDASetMaxStep(m_ida_mem, hmax);
156 checkError(flag, "setMaxStepSize", "IDASetMaxStep");
157 }
158}
159
161{
162 m_maxsteps = nmax;
163 if (m_ida_mem) {
164 IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
165 }
166}
167
169{
170 return m_maxsteps;
171}
172
174{
176 if (m_ida_mem) {
177 IDASetMaxErrTestFails(m_ida_mem, n);
178 }
179}
180
182{
183 AnyMap stats;
184 long int val;
185 int lastOrder;
186
187 int flag = IDAGetNumSteps(m_ida_mem, &val);
188 checkError(flag, "solverStats", "IDAGetNumSteps");
189 stats["steps"] = val;
190 IDAGetNumResEvals(m_ida_mem, &val);
191 stats["res_evals"] = val;
192 IDAGetNumLinSolvSetups(m_ida_mem, &val);
193 stats["lin_solve_setups"] = val;
194 IDAGetNumErrTestFails(m_ida_mem, &val);
195 stats["err_tests_fails"] = val;
196 IDAGetLastOrder(m_ida_mem, &lastOrder);
197 stats["last_order"] = lastOrder;
198 IDAGetNumNonlinSolvIters(m_ida_mem, &val);
199 stats["nonlinear_iters"] = val;
200 IDAGetNumNonlinSolvConvFails(m_ida_mem, &val);
201 stats["nonlinear_conv_fails"] = val;
202 return stats;
203}
204
205void IdasIntegrator::setMaxNonlinIterations(int n)
206{
208 if (m_ida_mem) {
209 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
210 checkError(flag, "setMaxNonlinIterations", "IDASetMaxNonlinIters");
211 }
212}
213
214void IdasIntegrator::setMaxNonlinConvFailures(int n)
215{
217 if (m_ida_mem) {
218 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
219 checkError(flag, "setMaxNonlinConvFailures", "IDASetMaxConvFails");
220 }
221}
222
223void IdasIntegrator::includeAlgebraicInErrorTest(bool yesno)
224{
225 m_setSuppressAlg = !yesno;
226 if (m_ida_mem) {
227 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
228 checkError(flag, "inclAlgebraicInErrorTest", "IDASetSuppressAlg");
229 }
230}
231
233{
234 m_neq = func.neq();
235 m_t0 = t0;
236 m_time = t0;
237 m_tInteg = t0;
238 m_func = &func;
239 func.clearErrors();
240
241 if (m_y) {
242 N_VDestroy_Serial(m_y); // free solution vector if already allocated
243 }
244 m_y = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
245 N_VConst(0.0, m_y);
246
247 if (m_ydot) {
248 N_VDestroy_Serial(m_ydot); // free derivative vector if already allocated
249 }
250 m_ydot = newNVector(m_neq, m_sundials_ctx);
251 N_VConst(0.0, m_ydot);
252
253 // check abs tolerance array size
254 if (m_itol == IDA_SV && m_nabs < m_neq) {
255 throw CanteraError("IdasIntegrator::initialize",
256 "not enough absolute tolerance values specified.");
257 }
258
259 if (m_constraints) {
260 N_VDestroy_Serial(m_constraints);
261 }
262 m_constraints = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
263 // set the constraints
264 func.getConstraints(NV_DATA_S(m_constraints));
265
266 // get the initial conditions
267 func.getStateDae(NV_DATA_S(m_y), NV_DATA_S(m_ydot));
268
269 if (m_ida_mem) {
270 IDAFree(&m_ida_mem);
271 }
272
273 //! Create the IDA solver
274 #if SUNDIALS_VERSION_MAJOR >= 6
275 m_ida_mem = IDACreate(m_sundials_ctx.get());
276 #else
277 m_ida_mem = IDACreate();
278 #endif
279 if (!m_ida_mem) {
280 throw CanteraError("IdasIntegrator::initialize", "IDACreate failed.");
281 }
282
283 int flag = IDAInit(m_ida_mem, ida_rhs, m_t0, m_y, m_ydot);
284 if (flag != IDA_SUCCESS) {
285 if (flag == IDA_MEM_FAIL) {
286 throw CanteraError("IdasIntegrator::initialize",
287 "Memory allocation failed.");
288 } else if (flag == IDA_ILL_INPUT) {
289 throw CanteraError("IdasIntegrator::initialize",
290 "Illegal value for IDAInit input argument.");
291 } else {
292 throw CanteraError("IdasIntegrator::initialize", "IDAInit failed.");
293 }
294 }
295
296 #if SUNDIALS_VERSION_MAJOR >= 7
297 flag = SUNContext_PushErrHandler(m_sundials_ctx.get(), &sundials_err, this);
298 #else
299 flag = IDASetErrHandlerFn(m_ida_mem, &ida_err, this);
300 #endif
301
302 // set constraints
303 flag = IDASetId(m_ida_mem, m_constraints);
304 checkError(flag, "initialize", "IDASetId");
305
306 if (m_itol == IDA_SV) {
307 flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
308 checkError(flag, "initialize", "IDASVtolerances");
309 } else {
310 flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
311 checkError(flag, "initialize", "IDASStolerances");
312 }
313
314 flag = IDASetUserData(m_ida_mem, &func);
315 checkError(flag, "initialize", "IDASetUserData");
316
317 if (func.nparams() > 0) {
318 throw CanteraError("IdasIntegrator::initialize", "Sensitivity analysis "
319 "for DAE systems is not fully implemented");
320 sensInit(t0, func);
321 flag = IDASetSensParams(m_ida_mem, func.m_sens_params.data(),
322 func.m_paramScales.data(), NULL);
323 checkError(flag, "initialize", "IDASetSensParams");
324 }
325 applyOptions();
326}
327
328void IdasIntegrator::reinitialize(double t0, FuncEval& func)
329{
330 m_t0 = t0;
331 m_time = t0;
332 m_tInteg = t0;
333 func.getStateDae(NV_DATA_S(m_y), NV_DATA_S(m_ydot));
334 m_func = &func;
335 func.clearErrors();
336
337 int result = IDAReInit(m_ida_mem, m_t0, m_y, m_ydot);
338 checkError(result, "reinitialize", "IDAReInit");
339 applyOptions();
340}
341
343{
344 if (m_type == "DENSE") {
345 sd_size_t N = static_cast<sd_size_t>(m_neq);
346 SUNLinSolFree((SUNLinearSolver) m_linsol);
347 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
348 #if SUNDIALS_VERSION_MAJOR >= 6
349 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
350 #else
351 m_linsol_matrix = SUNDenseMatrix(N, N);
352 #endif
353 #if SUNDIALS_VERSION_MAJOR >= 6
354 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
355 #else
356 m_linsol_matrix = SUNDenseMatrix(N, N);
357 #endif
358 #if CT_SUNDIALS_USE_LAPACK
359 #if SUNDIALS_VERSION_MAJOR >= 6
360 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
361 m_sundials_ctx.get());
362 #else
363 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
364 #endif
365 #else
366 #if SUNDIALS_VERSION_MAJOR >= 6
367 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
368 m_sundials_ctx.get());
369 #else
370 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix);
371 #endif
372 #endif
373 #if SUNDIALS_VERSION_MAJOR >= 6
374 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
375 (SUNMatrix) m_linsol_matrix);
376 #else
377 IDADlsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol,
378 (SUNMatrix) m_linsol_matrix);
379 #endif
380 } else if (m_type == "GMRES") {
381 SUNLinSolFree((SUNLinearSolver) m_linsol);
382 #if SUNDIALS_VERSION_MAJOR >= 6
383 m_linsol = SUNLinSol_SPGMR(m_y, SUN_PREC_NONE, 0, m_sundials_ctx.get());
384 IDASetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol, nullptr);
385 #elif SUNDIALS_VERSION_MAJOR >= 4
386 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
387 IDASpilsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol);
388 #else
389 m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
390 IDASpilsSetLinearSolver(m_ida_mem, (SUNLinearSolver) m_linsol);
391 #endif
392 } else {
393 throw CanteraError("IdasIntegrator::applyOptions",
394 "unsupported linear solver flag '{}'", m_type);
395 }
396
397 if (m_init_step > 0) {
398 IDASetInitStep(m_ida_mem, m_init_step);
399 }
400
401 if (m_maxord > 0) {
402 int flag = IDASetMaxOrd(m_ida_mem, m_maxord);
403 checkError(flag, "applyOptions", "IDASetMaxOrd");
404 }
405 if (m_maxsteps > 0) {
406 IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
407 }
408 if (m_hmax > 0) {
409 IDASetMaxStep(m_ida_mem, m_hmax);
410 }
411 if (m_maxNonlinIters > 0) {
412 int flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
413 checkError(flag, "applyOptions", "IDASetMaxNonlinIters");
414 }
415 if (m_maxNonlinConvFails > 0) {
416 int flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
417 checkError(flag, "applyOptions", "IDASetMaxConvFails");
418 }
419 if (m_setSuppressAlg != 0) {
420 int flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
421 checkError(flag, "applyOptions", "IDASetSuppressAlg");
422 }
423}
424
425void IdasIntegrator::sensInit(double t0, FuncEval& func)
426{
427 m_np = func.nparams();
428 m_sens_ok = false;
429
430 N_Vector y = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
431 #if SUNDIALS_VERSION_MAJOR >= 6
432 m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
433 #else
434 m_yS = N_VCloneVectorArray_Serial(static_cast<int>(m_np), y);
435 #endif
436 for (size_t n = 0; n < m_np; n++) {
437 N_VConst(0.0, m_yS[n]);
438 }
439 N_VDestroy_Serial(y);
440 N_Vector ydot = newNVector(static_cast<sd_size_t>(func.neq()), m_sundials_ctx);
441 #if SUNDIALS_VERSION_MAJOR >= 6
442 m_ySdot = N_VCloneVectorArray(static_cast<int>(m_np), ydot);
443 #else
444 m_ySdot = N_VCloneVectorArray_Serial(static_cast<int>(m_np), ydot);
445 #endif
446 for (size_t n = 0; n < m_np; n++) {
447 N_VConst(0.0, m_ySdot[n]);
448 }
449
450 int flag = IDASensInit(m_ida_mem, static_cast<sd_size_t>(m_np),
451 IDA_STAGGERED, IDASensResFn(0), m_yS, m_ySdot);
452 checkError(flag, "sensInit", "IDASensInit");
453
454 vector<double> atol(m_np);
455 for (size_t n = 0; n < m_np; n++) {
456 // This scaling factor is tuned so that reaction and species enthalpy
457 // sensitivities can be computed simultaneously with the same abstol.
458 atol[n] = m_abstolsens / func.m_paramScales[n];
459 }
460 flag = IDASensSStolerances(m_ida_mem, m_reltolsens, atol.data());
461 checkError(flag, "sensInit", "IDASensSStolerances");
462}
463
465{
466 if (tout == m_time) {
467 return;
468 } else if (tout < m_time) {
469 throw CanteraError("IdasIntegrator::integrate",
470 "Cannot integrate backwards in time.\n"
471 "Requested time {} < current time {}",
472 tout, m_time);
473 }
474 int nsteps = 0;
475 while (m_tInteg < tout) {
476 if (nsteps >= m_maxsteps) {
477 throw CanteraError("IdasIntegrator::integrate",
478 "Maximum number of timesteps ({}) taken without reaching output "
479 "time ({}).\nCurrent integrator time: {}",
480 nsteps, tout, m_time);
481 }
482 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
483 if (flag != IDA_SUCCESS) {
484 string f_errs = m_func->getErrors();
485 if (!f_errs.empty()) {
486 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
487 }
488 throw CanteraError("IdasIntegrator::integrate",
489 "IDA error encountered. Error code: {}\n{}\n"
490 "{}"
491 "Components with largest weighted error estimates:\n{}",
492 flag, m_error_message, f_errs, getErrorInfo(10));
493 }
494 nsteps++;
495 }
496 int flag = IDAGetDky(m_ida_mem, tout, 0, m_y);
497 checkError(flag, "integrate", "IDAGetDky");
498 m_time = tout;
499}
500
501double IdasIntegrator::step(double tout)
502{
503 int flag = IDASolve(m_ida_mem, tout, &m_tInteg, m_y, m_ydot, IDA_ONE_STEP);
504 if (flag != IDA_SUCCESS) {
505 string f_errs = m_func->getErrors();
506 if (!f_errs.empty()) {
507 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
508 }
509 throw CanteraError("IdasIntegrator::step",
510 "IDA error encountered. Error code: {}\n{}\n"
511 "{}"
512 "Components with largest weighted error estimates:\n{}",
513 flag, f_errs, m_error_message, getErrorInfo(10));
514
515 }
517 return m_time;
518}
519
520double IdasIntegrator::sensitivity(size_t k, size_t p)
521{
522 if (m_time == m_t0) {
523 // calls to IDAGetSens are only allowed after a successful time step.
524 return 0.0;
525 }
526 if (!m_sens_ok && m_np) {
527 int flag = IDAGetSensDky(m_ida_mem, m_time, 0, m_yS);
528 checkError(flag, "sensitivity", "IDAGetSens");
529 m_sens_ok = true;
530 }
531
532 if (k >= m_neq) {
533 throw CanteraError("IdasIntegrator::sensitivity",
534 "sensitivity: k out of range ({})", k);
535 }
536 if (p >= m_np) {
537 throw CanteraError("IdasIntegrator::sensitivity",
538 "sensitivity: p out of range ({})", p);
539 }
540 return NV_Ith_S(m_yS[p],k);
541}
542
544{
545 N_Vector errs = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
546 N_Vector errw = newNVector(static_cast<sd_size_t>(m_neq), m_sundials_ctx);
547 IDAGetErrWeights(m_ida_mem, errw);
548 IDAGetEstLocalErrors(m_ida_mem, errs);
549
550 vector<tuple<double, double, size_t>> weightedErrors;
551 for (size_t i=0; i<m_neq; i++) {
552 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
553 weightedErrors.emplace_back(-abs(err), err, i);
554 }
555 N_VDestroy(errs);
556 N_VDestroy(errw);
557
558 N = std::min(N, static_cast<int>(m_neq));
559 sort(weightedErrors.begin(), weightedErrors.end());
560 fmt::memory_buffer s;
561 for (int i=0; i<N; i++) {
562 fmt_append(s, "{}: {}\n", get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
563 }
564 return to_string(s);
565}
566
567void IdasIntegrator::checkError(long flag, const string& ctMethod,
568 const string& idaMethod) const
569{
570 if (flag == IDA_SUCCESS) {
571 return;
572 } else if (flag == IDA_MEM_NULL) {
573 throw CanteraError("IdasIntegrator::" + ctMethod,
574 "IDAS integrator is not initialized");
575 } else {
576 const char* flagname = IDAGetReturnFlagName(flag);
577 throw CanteraError("IdasIntegrator::" + ctMethod,
578 "{} returned error code {} ({}):\n{}",
579 idaMethod, flag, flagname, m_error_message);
580 }
581}
582
584{
585 if (t != BDF_Method) {
586 // IDA only has the BDF method
587 throw CanteraError("IdasIntegrator::setMethod", "unknown method");
588 }
589}
590
591}
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(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, const std::string &tmpl, 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: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
Contains declarations for string manipulation functions within Cantera.