Cantera  3.2.0
Loading...
Searching...
No Matches
CVodesIntegrator.cpp
Go to the documentation of this file.
1//! @file CVodesIntegrator.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 <iostream>
10using namespace std;
11
12#include "cantera/numerics/sundials_headers.h"
13
14namespace {
15
16N_Vector newNVector(size_t N, Cantera::SundialsContext& context)
17{
18#if SUNDIALS_VERSION_MAJOR >= 6
19 return N_VNew_Serial(static_cast<sd_size_t>(N), context.get());
20#else
21 return N_VNew_Serial(static_cast<sd_size_t>(N));
22#endif
23}
24
25} // end anonymous namespace
26
27namespace Cantera
28{
29
30extern "C" {
31 /**
32 * Function called by cvodes to evaluate ydot given y. The CVODE integrator
33 * allows passing in a void* pointer to access external data. This pointer
34 * is cast to a pointer to a instance of class FuncEval. The equations to be
35 * integrated should be specified by deriving a class from FuncEval that
36 * evaluates the desired equations.
37 * @ingroup odeGroup
38 */
39 static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void* f_data)
40 {
41 FuncEval* f = (FuncEval*) f_data;
42 return f->evalNoThrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
43 }
44
45 #if SUNDIALS_VERSION_MAJOR >= 7
46 //! Function called by CVodes when an error is encountered instead of
47 //! writing to stdout. Here, save the error message provided by CVodes so
48 //! that it can be included in the subsequently raised CanteraError. Used by
49 //! SUNDIALS 7.0 and newer.
50 static void sundials_err(int line, const char *func, const char *file,
51 const char *msg, SUNErrCode err_code,
52 void *err_user_data, SUNContext sunctx)
53 {
54 CVodesIntegrator* integrator = (CVodesIntegrator*) err_user_data;
55 integrator->m_error_message = fmt::format("{}: {}\n", func, msg);
56 }
57 #else
58 //! Function called by CVodes when an error is encountered instead of
59 //! writing to stdout. Here, save the error message provided by CVodes so
60 //! that it can be included in the subsequently raised CanteraError. Used by
61 //! SUNDIALS 6.x and older.
62 static void cvodes_err(int error_code, const char* module,
63 const char* function, char* msg, void* eh_data)
64 {
65 CVodesIntegrator* integrator = (CVodesIntegrator*) eh_data;
66 integrator->m_error_message = msg;
67 integrator->m_error_message += "\n";
68 }
69 #endif
70
71 static int cvodes_prec_setup(sunrealtype t, N_Vector y, N_Vector ydot,
72 sunbooleantype jok, sunbooleantype *jcurPtr,
73 sunrealtype gamma, void *f_data)
74 {
75 FuncEval* f = (FuncEval*) f_data;
76 if (!jok) {
77 *jcurPtr = true; // jacobian data was recomputed
78 return f->preconditioner_setup_nothrow(t, NV_DATA_S(y), gamma);
79 } else {
80 f->updatePreconditioner(gamma); // updates preconditioner with new gamma
81 *jcurPtr = false; // indicates that Jacobian data was not recomputed
82 return 0; // no error because not recomputed
83 }
84 }
85
86 static int cvodes_prec_solve(sunrealtype t, N_Vector y, N_Vector ydot, N_Vector r,
87 N_Vector z, sunrealtype gamma, sunrealtype delta,
88 int lr, void* f_data)
89 {
90 FuncEval* f = (FuncEval*) f_data;
91 return f->preconditioner_solve_nothrow(NV_DATA_S(r),NV_DATA_S(z));
92 }
93
94 /**
95 * SUNDIALS callback that forwards root evaluations to the FuncEval.
96 * @param[in] t Current integration time at which roots are requested
97 * @param[in] y Solution vector provided by CVODE (length FuncEval::neq())
98 * @param[out] gout Output array that must be filled with FuncEval::nRootFunctions()
99 * values
100 * @param[in] user_data Opaque pointer (FuncEval*) supplied during integrator setup
101 * @returns The result of FuncEval::evalRootFunctions (0 on success)
102 */
103 static int cvodes_root(sunrealtype t, N_Vector y, sunrealtype* gout, void* user_data)
104 {
105 auto* f = static_cast<FuncEval*>(user_data);
106 return f->evalRootFunctionsNoThrow(t, NV_DATA_S(y), gout);
107 }
108}
109
111 : m_itol(CV_SS)
112 , m_method(CV_BDF)
113{
114}
115
116CVodesIntegrator::~CVodesIntegrator()
117{
118 if (m_cvode_mem) {
119 if (m_np > 0) {
120 CVodeSensFree(m_cvode_mem);
121 }
122 CVodeFree(&m_cvode_mem);
123 }
124
125 SUNLinSolFree((SUNLinearSolver) m_linsol);
126 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
127
128 if (m_y) {
129 N_VDestroy_Serial(m_y);
130 }
131 if (m_abstol) {
132 N_VDestroy_Serial(m_abstol);
133 }
134 if (m_dky) {
135 N_VDestroy_Serial(m_dky);
136 }
137 if (m_yS) {
138 #if SUNDIALS_VERSION_MAJOR >= 6
139 N_VDestroyVectorArray(m_yS, static_cast<int>(m_np));
140 #else
141 N_VDestroyVectorArray_Serial(m_yS, static_cast<int>(m_np));
142 #endif
143 }
144}
145
146
148{
149 return NV_Ith_S(m_y, k);
150}
151
153{
154 return NV_DATA_S(m_y);
155}
156
157void CVodesIntegrator::setTolerances(double reltol, size_t n, double* abstol)
158{
159 m_itol = CV_SV;
160 m_nabs = n;
161 if (n != m_neq) {
162 if (m_abstol) {
163 N_VDestroy_Serial(m_abstol);
164 }
165 m_abstol = newNVector(n, m_sundials_ctx);
166 }
167 for (size_t i=0; i<n; i++) {
168 NV_Ith_S(m_abstol, i) = abstol[i];
169 }
170 m_reltol = reltol;
171}
172
173void CVodesIntegrator::setTolerances(double reltol, double abstol)
174{
175 m_itol = CV_SS;
176 m_reltol = reltol;
177 m_abstols = abstol;
178}
179
180void CVodesIntegrator::setSensitivityTolerances(double reltol, double abstol)
181{
182 m_reltolsens = reltol;
183 m_abstolsens = abstol;
184}
185
187{
188 if (t == BDF_Method) {
189 m_method = CV_BDF;
190 } else if (t == Adams_Method) {
191 m_method = CV_ADAMS;
192 } else {
193 throw CanteraError("CVodesIntegrator::setMethod", "unknown method");
194 }
195}
196
198{
199 m_hmax = hmax;
200 if (m_cvode_mem) {
201 CVodeSetMaxStep(m_cvode_mem, hmax);
202 }
203}
204
206{
207 m_hmin = hmin;
208 if (m_cvode_mem) {
209 CVodeSetMinStep(m_cvode_mem, hmin);
210 }
211}
212
214{
215 m_maxsteps = nmax;
216 if (m_cvode_mem) {
217 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
218 }
219}
220
222{
223 return m_maxsteps;
224}
225
227{
228 m_maxErrTestFails = n;
229 if (m_cvode_mem) {
230 CVodeSetMaxErrTestFails(m_cvode_mem, n);
231 }
232}
233
235{
236 // Skip the Sundials call when the requested count matches what CVODE already has
237 if (m_cvode_mem && m_nRootFunctions == nroots) {
238 return;
239 }
240 m_nRootFunctions = nroots;
241 // When initialize() hasn’t created CVODE memory yet, just cache the count
242 if (!m_cvode_mem) {
243 return;
244 }
245 // Register (or remove) the root callback; passing nullptr disables root finding
246 CVRootFn root_cb = nroots ? &cvodes_root : nullptr;
247 int flag = CVodeRootInit(m_cvode_mem, static_cast<int>(nroots), root_cb);
248 checkError(flag, "setRootFunctionCount", "CVodeRootInit");
249}
250
251void CVodesIntegrator::sensInit(double t0, FuncEval& func)
252{
253 m_np = func.nparams();
254 m_sens_ok = false;
255
256 N_Vector y = newNVector(func.neq(), m_sundials_ctx);
257 #if SUNDIALS_VERSION_MAJOR >= 6
258 m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
259 #else
260 m_yS = N_VCloneVectorArray_Serial(static_cast<int>(m_np), y);
261 #endif
262 for (size_t n = 0; n < m_np; n++) {
263 N_VConst(0.0, m_yS[n]);
264 }
265 N_VDestroy_Serial(y);
266
267 int flag = CVodeSensInit(m_cvode_mem, static_cast<int>(m_np),
268 CV_STAGGERED, CVSensRhsFn(0), m_yS);
269 checkError(flag, "sensInit", "CVodeSensInit");
270
271 vector<double> atol(m_np);
272 for (size_t n = 0; n < m_np; n++) {
273 // This scaling factor is tuned so that reaction and species enthalpy
274 // sensitivities can be computed simultaneously with the same abstol.
275 atol[n] = m_abstolsens / func.m_paramScales[n];
276 }
277 flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
278 checkError(flag, "sensInit", "CVodeSensSStolerances");
279}
280
282{
283 m_neq = func.neq();
284 m_t0 = t0;
285 m_time = t0;
286 m_tInteg = t0;
287 m_func = &func;
288 func.clearErrors();
289 // Initialize preconditioner if applied
290 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
291 m_preconditioner->initialize(m_neq);
292 }
293 if (m_y) {
294 N_VDestroy_Serial(m_y); // free solution vector if already allocated
295 }
296 m_y = newNVector(m_neq, m_sundials_ctx); // allocate solution vector
297 N_VConst(0.0, m_y);
298 if (m_dky) {
299 N_VDestroy_Serial(m_dky); // free derivative vector if already allocated
300 }
301 m_dky = newNVector(m_neq, m_sundials_ctx); // allocate derivative vector
302 N_VConst(0.0, m_dky);
303 // check abs tolerance array size
304 if (m_itol == CV_SV && m_nabs < m_neq) {
305 throw CanteraError("CVodesIntegrator::initialize",
306 "not enough absolute tolerance values specified.");
307 }
308
309 func.getState(NV_DATA_S(m_y));
310
311 if (m_cvode_mem) {
312 CVodeFree(&m_cvode_mem);
313 }
314
315 //! Specify the method and the iteration type. Cantera Defaults:
316 //! CV_BDF - Use BDF methods
317 //! CV_NEWTON - use Newton's method
318 #if SUNDIALS_VERSION_MAJOR < 6
319 m_cvode_mem = CVodeCreate(m_method);
320 #else
321 m_cvode_mem = CVodeCreate(m_method, m_sundials_ctx.get());
322 #endif
323 if (!m_cvode_mem) {
324 throw CanteraError("CVodesIntegrator::initialize",
325 "CVodeCreate failed.");
326 }
327
328 int flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, m_y);
329 if (flag != CV_SUCCESS) {
330 if (flag == CV_MEM_FAIL) {
331 throw CanteraError("CVodesIntegrator::initialize",
332 "Memory allocation failed.");
333 } else if (flag == CV_ILL_INPUT) {
334 throw CanteraError("CVodesIntegrator::initialize",
335 "Illegal value for CVodeInit input argument.");
336 } else {
337 throw CanteraError("CVodesIntegrator::initialize",
338 "CVodeInit failed.");
339 }
340 }
341 #if SUNDIALS_VERSION_MAJOR >= 7
342 flag = SUNContext_PushErrHandler(m_sundials_ctx.get(), &sundials_err, this);
343 #else
344 flag = CVodeSetErrHandlerFn(m_cvode_mem, &cvodes_err, this);
345 #endif
346
347 if (m_itol == CV_SV) {
348 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
349 checkError(flag, "initialize", "CVodeSVtolerances");
350 } else {
351 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
352 checkError(flag, "initialize", "CVodeSStolerances");
353 }
354
355 flag = CVodeSetUserData(m_cvode_mem, &func);
356 checkError(flag, "initialize", "CVodeSetUserData");
357
358 m_nRootFunctions = npos;
360
361 if (func.nparams() > 0) {
362 sensInit(t0, func);
363 flag = CVodeSetSensParams(m_cvode_mem, func.m_sens_params.data(),
364 func.m_paramScales.data(), NULL);
365 checkError(flag, "initialize", "CVodeSetSensParams");
366 }
367 applyOptions();
368}
369
370void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
371{
372 m_t0 = t0;
373 m_time = t0;
374 m_tInteg = t0;
375 func.getState(NV_DATA_S(m_y));
376 m_func = &func;
377 func.clearErrors();
378 // reinitialize preconditioner if applied
379 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
380 m_preconditioner->initialize(m_neq);
381 }
382 int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
383 checkError(result, "reinitialize", "CVodeReInit");
384 m_nRootFunctions = npos;
386 applyOptions();
387}
388
390{
391 if (m_type == "DENSE") {
392 sd_size_t N = static_cast<sd_size_t>(m_neq);
393 SUNLinSolFree((SUNLinearSolver) m_linsol);
394 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
395 #if SUNDIALS_VERSION_MAJOR >= 6
396 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
397 #else
398 m_linsol_matrix = SUNDenseMatrix(N, N);
399 #endif
400 if (m_linsol_matrix == nullptr) {
401 throw CanteraError("CVodesIntegrator::applyOptions",
402 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
403 }
404 int flag;
405 #if SUNDIALS_VERSION_MAJOR >= 6
406 #if CT_SUNDIALS_USE_LAPACK
407 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
408 m_sundials_ctx.get());
409 #else
410 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
411 m_sundials_ctx.get());
412 #endif
413 flag = CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
414 (SUNMatrix) m_linsol_matrix);
415 #else
416 #if CT_SUNDIALS_USE_LAPACK
417 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
418 #else
419 m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
420 #endif
421 flag = CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
422 (SUNMatrix) m_linsol_matrix);
423 #endif
424 if (m_linsol == nullptr) {
425 throw CanteraError("CVodesIntegrator::applyOptions",
426 "Error creating Sundials dense linear solver object");
427 } else if (flag != CV_SUCCESS) {
428 throw CanteraError("CVodesIntegrator::applyOptions",
429 "Error connecting linear solver to CVODES. "
430 "Sundials error code: {}", flag);
431 }
432
433 // throw preconditioner error for DENSE + NOJAC
434 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
435 throw CanteraError("CVodesIntegrator::applyOptions",
436 "Preconditioning is not available with the specified problem type.");
437 }
438 } else if (m_type == "DIAG") {
439 CVDiag(m_cvode_mem);
440 // throw preconditioner error for DIAG
441 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
442 throw CanteraError("CVodesIntegrator::applyOptions",
443 "Preconditioning is not available with the specified problem type.");
444 }
445 } else if (m_type == "GMRES") {
446 SUNLinSolFree((SUNLinearSolver) m_linsol);
447 #if SUNDIALS_VERSION_MAJOR >= 6
448 m_linsol = SUNLinSol_SPGMR(m_y, SUN_PREC_NONE, 0, m_sundials_ctx.get());
449 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol, nullptr);
450 #else
451 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
452 CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
453 #endif
454 // set preconditioner if used
455 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
456 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver) m_linsol,
457 static_cast<int>(m_prec_side));
458 CVodeSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
459 cvodes_prec_solve);
460 }
461 } else if (m_type == "BAND") {
462 sd_size_t N = static_cast<sd_size_t>(m_neq);
463 sd_size_t nu = m_mupper;
464 sd_size_t nl = m_mlower;
465 SUNLinSolFree((SUNLinearSolver) m_linsol);
466 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
467 #if SUNDIALS_VERSION_MAJOR >= 6
468 m_linsol_matrix = SUNBandMatrix(N, nu, nl, m_sundials_ctx.get());
469 #else
470 m_linsol_matrix = SUNBandMatrix(N, nu, nl);
471 #endif
472 if (m_linsol_matrix == nullptr) {
473 throw CanteraError("CVodesIntegrator::applyOptions",
474 "Unable to create SUNBandMatrix of size {} with bandwidths "
475 "{} and {}", N, nu, nl);
476 }
477 #if SUNDIALS_VERSION_MAJOR >= 6
478 #if CT_SUNDIALS_USE_LAPACK
479 m_linsol = SUNLinSol_LapackBand(m_y, (SUNMatrix) m_linsol_matrix,
480 m_sundials_ctx.get());
481 #else
482 m_linsol = SUNLinSol_Band(m_y, (SUNMatrix) m_linsol_matrix,
483 m_sundials_ctx.get());
484 #endif
485 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
486 (SUNMatrix) m_linsol_matrix);
487 #else
488 #if CT_SUNDIALS_USE_LAPACK
489 m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
490 #else
491 m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
492 #endif
493 CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
494 (SUNMatrix) m_linsol_matrix);
495 #endif
496 } else {
497 throw CanteraError("CVodesIntegrator::applyOptions",
498 "unsupported linear solver flag '{}'", m_type);
499 }
500
501 if (m_maxord > 0) {
502 int flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
503 checkError(flag, "applyOptions", "CVodeSetMaxOrd");
504 }
505 if (m_maxsteps > 0) {
506 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
507 }
508 if (m_hmax > 0) {
509 CVodeSetMaxStep(m_cvode_mem, m_hmax);
510 }
511 if (m_hmin > 0) {
512 CVodeSetMinStep(m_cvode_mem, m_hmin);
513 }
514 if (m_maxErrTestFails > 0) {
515 CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
516 }
517}
518
520{
521 if (tout == m_time) {
522 return;
523 } else if (tout < m_time) {
524 throw CanteraError("CVodesIntegrator::integrate",
525 "Cannot integrate backwards in time.\n"
526 "Requested time {} < current time {}",
527 tout, m_time);
528 }
529 int nsteps = 0;
530 while (m_tInteg < tout) {
531 if (nsteps >= m_maxsteps) {
532 string f_errs = m_func->getErrors();
533 if (!f_errs.empty()) {
534 f_errs = "\nExceptions caught during RHS evaluation:\n" + f_errs;
535 }
536 throw CanteraError("CVodesIntegrator::integrate",
537 "Maximum number of timesteps ({}) taken without reaching output "
538 "time ({}).\nCurrent integrator time: {}{}",
539 nsteps, tout, m_tInteg, f_errs);
540 }
541 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
542 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
543 string f_errs = m_func->getErrors();
544 if (!f_errs.empty()) {
545 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
546 }
547 throw CanteraError("CVodesIntegrator::integrate",
548 "CVodes error encountered. Error code: {}\n{}\n"
549 "{}"
550 "Components with largest weighted error estimates:\n{}",
551 flag, m_error_message, f_errs, getErrorInfo(10));
552 }
553 if (flag == CV_ROOT_RETURN) {
554 // Stop early at root (e.g., advance limit reached); align tout to the
555 // root time
556 tout = m_tInteg;
557 break;
558 }
559 nsteps++;
560 }
561
562 // Interpolate the solution to either the user-specified output time or
563 // the time at which a root event occurred.
564 double t_eval = tout;
565 int flag = CVodeGetDky(m_cvode_mem, t_eval, 0, m_y);
566 checkError(flag, "integrate", "CVodeGetDky");
567 m_time = t_eval;
568 m_sens_ok = false;
569}
570
571double CVodesIntegrator::step(double tout)
572{
573 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
574 if (flag != CV_SUCCESS && flag != CV_ROOT_RETURN) {
575 string f_errs = m_func->getErrors();
576 if (!f_errs.empty()) {
577 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
578 }
579 throw CanteraError("CVodesIntegrator::step",
580 "CVodes error encountered. Error code: {}\n{}\n"
581 "{}"
582 "Components with largest weighted error estimates:\n{}",
583 flag, f_errs, m_error_message, getErrorInfo(10));
584
585 }
586 m_sens_ok = false;
588 return m_time;
589}
590
591double* CVodesIntegrator::derivative(double tout, int n)
592{
593 int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
594 checkError(flag, "derivative", "CVodeGetDky");
595 return NV_DATA_S(m_dky);
596}
597
599{
600 int ord;
601 CVodeGetLastOrder(m_cvode_mem, &ord);
602 return ord;
603}
604
606{
607 long int ne;
608 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
609 return ne;
610}
611
613{
614 // AnyMap to return stats
615 AnyMap stats;
616
617 // long int linear solver stats provided by CVodes
618 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
619 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
620 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
621 nonlinConvFails = 0, orderReductions = 0;
622 int lastOrder = 0;
623;
624
625 CVodeGetNumSteps(m_cvode_mem, &steps);
626 CVodeGetNumRhsEvals(m_cvode_mem, &rhsEvals);
627 CVodeGetNonlinSolvStats(m_cvode_mem, &nonlinIters, &nonlinConvFails);
628 CVodeGetNumErrTestFails(m_cvode_mem, &errTestFails);
629 CVodeGetLastOrder(m_cvode_mem, &lastOrder);
630 CVodeGetNumStabLimOrderReds(m_cvode_mem, &orderReductions);
631 CVodeGetNumJacEvals(m_cvode_mem, &jacEvals);
632 CVodeGetNumLinRhsEvals(m_cvode_mem, &linRhsEvals);
633 CVodeGetNumLinSolvSetups(m_cvode_mem, &linSetup);
634 CVodeGetNumLinIters(m_cvode_mem, &linIters);
635 CVodeGetNumLinConvFails(m_cvode_mem, &linConvFails);
636 CVodeGetNumPrecEvals(m_cvode_mem, &precEvals);
637 CVodeGetNumPrecSolves(m_cvode_mem, &precSolves);
638 CVodeGetNumJTSetupEvals(m_cvode_mem, &jtSetupEvals);
639 CVodeGetNumJtimesEvals(m_cvode_mem, &jTimesEvals);
640
641 #if SUNDIALS_VERSION_MAJOR >= 7 || (SUNDIALS_VERSION_MAJOR == 6 && SUNDIALS_VERSION_MINOR >= 2)
642 long int stepSolveFails = 0;
643 CVodeGetNumStepSolveFails(m_cvode_mem, &stepSolveFails);
644 stats["step_solve_fails"] = stepSolveFails;
645 #endif
646
647 stats["steps"] = steps;
648 stats["rhs_evals"] = rhsEvals;
649 stats["nonlinear_iters"] = nonlinIters;
650 stats["nonlinear_conv_fails"] = nonlinConvFails;
651 stats["err_test_fails"] = errTestFails;
652 stats["last_order"] = lastOrder;
653 stats["stab_order_reductions"] = orderReductions;
654
655 stats["jac_evals"] = jacEvals;
656 stats["lin_solve_setups"] = linSetup;
657 stats["lin_rhs_evals"] = linRhsEvals;
658 stats["lin_iters"] = linIters;
659 stats["lin_conv_fails"] = linConvFails;
660 stats["prec_evals"] = precEvals;
661 stats["prec_solves"] = precSolves;
662 stats["jt_vec_setup_evals"] = jtSetupEvals;
663 stats["jt_vec_prod_evals"] = jTimesEvals;
664 return stats;
665}
666
667double CVodesIntegrator::sensitivity(size_t k, size_t p)
668{
669 if (m_time == m_t0) {
670 // calls to CVodeGetSens are only allowed after a successful time step.
671 return 0.0;
672 }
673 if (!m_sens_ok && m_np) {
674 int flag = CVodeGetSensDky(m_cvode_mem, m_time, 0, m_yS);
675 checkError(flag, "sensitivity", "CVodeGetSens");
676 m_sens_ok = true;
677 }
678
679 if (k >= m_neq) {
680 throw CanteraError("CVodesIntegrator::sensitivity",
681 "sensitivity: k out of range ({})", k);
682 }
683 if (p >= m_np) {
684 throw CanteraError("CVodesIntegrator::sensitivity",
685 "sensitivity: p out of range ({})", p);
686 }
687 return NV_Ith_S(m_yS[p],k);
688}
689
691{
692 N_Vector errs = newNVector(m_neq, m_sundials_ctx);
693 N_Vector errw = newNVector(m_neq, m_sundials_ctx);
694 CVodeGetErrWeights(m_cvode_mem, errw);
695 CVodeGetEstLocalErrors(m_cvode_mem, errs);
696
697 vector<tuple<double, double, size_t>> weightedErrors;
698 for (size_t i=0; i<m_neq; i++) {
699 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
700 weightedErrors.emplace_back(-abs(err), err, i);
701 }
702 N_VDestroy(errs);
703 N_VDestroy(errw);
704
705 N = std::min(N, static_cast<int>(m_neq));
706 sort(weightedErrors.begin(), weightedErrors.end());
707 fmt::memory_buffer s;
708 for (int i=0; i<N; i++) {
709 fmt_append(s, "{}: {}\n",
710 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
711 }
712 return to_string(s);
713}
714
715void CVodesIntegrator::checkError(long flag, const string& ctMethod,
716 const string& cvodesMethod) const
717{
718 if (flag == CV_SUCCESS) {
719 return;
720 } else if (flag == CV_MEM_NULL) {
721 throw CanteraError("CVodesIntegrator::" + ctMethod,
722 "CVODES integrator is not initialized");
723 } else {
724 const char* flagname = CVodeGetReturnFlagName(flag);
725 throw CanteraError("CVodesIntegrator::" + ctMethod,
726 "{} returned error code {} ({}):\n{}",
727 cvodesMethod, flag, flagname, m_error_message);
728 }
729}
730
731}
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
Wrapper class for 'cvodes' integrator from LLNL.
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.
N_Vector m_y
The system state at m_time.
void * m_cvode_mem
CVODES internal memory object.
void setRootFunctionCount(size_t nroots) override
Configure how many event/root functions the integrator should monitor.
void * m_linsol
Sundials linear solver object.
int nEvals() const override
The number of function evaluations.
double m_time
The current system time, corresponding to m_y.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS have been updated for at the current integrator ti...
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.
void checkError(long flag, const string &ctMethod, const string &cvodesMethod) const
Check whether a CVODES method indicated an error.
void setSensitivityTolerances(double reltol, double abstol) override
Set the sensitivity error tolerances.
void applyOptions()
Applies user-specified options to the underlying CVODES solver.
int lastOrder() const override
Order used during the last solution step.
void setMinStepSize(double hmin) override
Set the minimum step size.
void integrate(double tout) override
Integrate the system of equations.
void setTolerances(double reltol, size_t n, double *abstol) override
Set error tolerances.
double * derivative(double tout, int n) override
n-th derivative of the output function at time tout.
string m_error_message
Error message information provide by CVodes.
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.
AnyMap solverStats() const override
Get solver stats from integrator.
void * m_linsol_matrix
matrix used by Sundials
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.
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.
Base class for exceptions thrown by Cantera classes.
Virtual base class for ODE/DAE right-hand-side function evaluators.
Definition FuncEval.h:32
int evalRootFunctionsNoThrow(double t, const double *y, double *gout)
Wrapper for evalRootFunctions that converts exceptions to return codes.
Definition FuncEval.cpp:146
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
Definition FuncEval.h:206
void clearErrors()
Clear any previously-stored suppressed errors.
Definition FuncEval.h:196
int evalNoThrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
Definition FuncEval.cpp:7
virtual size_t neq() const =0
Number of equations.
virtual size_t nRootFunctions() const
Number of event/root functions exposed to the integrator.
Definition FuncEval.h:142
virtual size_t nparams() const
Number of sensitivity parameters.
Definition FuncEval.h:178
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:203
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
Definition FuncEval.h:164
shared_ptr< SystemJacobian > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
Definition Integrator.h:310
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
Definition Integrator.h:312
A wrapper for managing a SUNContext object, need for Sundials >= 6.0.
static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
MethodType
Specifies the method used to integrate the system of equations.
Definition Integrator.h:23
@ BDF_Method
Backward Differentiation.
Definition Integrator.h:24
@ Adams_Method
Adams.
Definition Integrator.h:25
static int cvodes_root(sunrealtype t, N_Vector y, sunrealtype *gout, void *user_data)
SUNDIALS callback that forwards root evaluations to the FuncEval.
static void cvodes_err(int error_code, const char *module, const char *function, char *msg, void *eh_data)
Function called by CVodes when an error is encountered instead of writing to stdout.
Contains declarations for string manipulation functions within Cantera.