Cantera  3.2.0a1
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 //! Function called by CVodes when an error is encountered instead of
46 //! writing to stdout. Here, save the error message provided by CVodes so
47 //! that it can be included in the subsequently raised CanteraError. Used by
48 //! SUNDIALS 6.x and older.
49 static void cvodes_err(int error_code, const char* module,
50 const char* function, char* msg, void* eh_data)
51 {
52 CVodesIntegrator* integrator = (CVodesIntegrator*) 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 CVodesIntegrator* integrator = (CVodesIntegrator*) err_user_data;
67 integrator->m_error_message = fmt::format("{}: {}\n", func, msg);
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
96 : m_itol(CV_SS)
97 , m_method(CV_BDF)
98{
99}
100
101CVodesIntegrator::~CVodesIntegrator()
102{
103 if (m_cvode_mem) {
104 if (m_np > 0) {
105 CVodeSensFree(m_cvode_mem);
106 }
107 CVodeFree(&m_cvode_mem);
108 }
109
110 SUNLinSolFree((SUNLinearSolver) m_linsol);
111 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
112
113 if (m_y) {
114 N_VDestroy_Serial(m_y);
115 }
116 if (m_abstol) {
117 N_VDestroy_Serial(m_abstol);
118 }
119 if (m_dky) {
120 N_VDestroy_Serial(m_dky);
121 }
122 if (m_yS) {
123 #if SUNDIALS_VERSION_MAJOR >= 6
124 N_VDestroyVectorArray(m_yS, static_cast<int>(m_np));
125 #else
126 N_VDestroyVectorArray_Serial(m_yS, static_cast<int>(m_np));
127 #endif
128 }
129}
130
131
133{
134 return NV_Ith_S(m_y, k);
135}
136
138{
139 return NV_DATA_S(m_y);
140}
141
142void CVodesIntegrator::setTolerances(double reltol, size_t n, double* abstol)
143{
144 m_itol = CV_SV;
145 m_nabs = n;
146 if (n != m_neq) {
147 if (m_abstol) {
148 N_VDestroy_Serial(m_abstol);
149 }
150 m_abstol = newNVector(n, m_sundials_ctx);
151 }
152 for (size_t i=0; i<n; i++) {
153 NV_Ith_S(m_abstol, i) = abstol[i];
154 }
155 m_reltol = reltol;
156}
157
158void CVodesIntegrator::setTolerances(double reltol, double abstol)
159{
160 m_itol = CV_SS;
161 m_reltol = reltol;
162 m_abstols = abstol;
163}
164
165void CVodesIntegrator::setSensitivityTolerances(double reltol, double abstol)
166{
167 m_reltolsens = reltol;
168 m_abstolsens = abstol;
169}
170
172{
173 if (t == BDF_Method) {
174 m_method = CV_BDF;
175 } else if (t == Adams_Method) {
176 m_method = CV_ADAMS;
177 } else {
178 throw CanteraError("CVodesIntegrator::setMethod", "unknown method");
179 }
180}
181
183{
184 m_hmax = hmax;
185 if (m_cvode_mem) {
186 CVodeSetMaxStep(m_cvode_mem, hmax);
187 }
188}
189
191{
192 m_hmin = hmin;
193 if (m_cvode_mem) {
194 CVodeSetMinStep(m_cvode_mem, hmin);
195 }
196}
197
199{
200 m_maxsteps = nmax;
201 if (m_cvode_mem) {
202 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
203 }
204}
205
207{
208 return m_maxsteps;
209}
210
212{
213 m_maxErrTestFails = n;
214 if (m_cvode_mem) {
215 CVodeSetMaxErrTestFails(m_cvode_mem, n);
216 }
217}
218
219void CVodesIntegrator::sensInit(double t0, FuncEval& func)
220{
221 m_np = func.nparams();
222 m_sens_ok = false;
223
224 N_Vector y = newNVector(func.neq(), m_sundials_ctx);
225 #if SUNDIALS_VERSION_MAJOR >= 6
226 m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
227 #else
228 m_yS = N_VCloneVectorArray_Serial(static_cast<int>(m_np), y);
229 #endif
230 for (size_t n = 0; n < m_np; n++) {
231 N_VConst(0.0, m_yS[n]);
232 }
233 N_VDestroy_Serial(y);
234
235 int flag = CVodeSensInit(m_cvode_mem, static_cast<int>(m_np),
236 CV_STAGGERED, CVSensRhsFn(0), m_yS);
237 checkError(flag, "sensInit", "CVodeSensInit");
238
239 vector<double> atol(m_np);
240 for (size_t n = 0; n < m_np; n++) {
241 // This scaling factor is tuned so that reaction and species enthalpy
242 // sensitivities can be computed simultaneously with the same abstol.
243 atol[n] = m_abstolsens / func.m_paramScales[n];
244 }
245 flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
246 checkError(flag, "sensInit", "CVodeSensSStolerances");
247}
248
250{
251 m_neq = func.neq();
252 m_t0 = t0;
253 m_time = t0;
254 m_tInteg = t0;
255 m_func = &func;
256 func.clearErrors();
257 // Initialize preconditioner if applied
258 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
259 m_preconditioner->initialize(m_neq);
260 }
261 if (m_y) {
262 N_VDestroy_Serial(m_y); // free solution vector if already allocated
263 }
264 m_y = newNVector(m_neq, m_sundials_ctx); // allocate solution vector
265 N_VConst(0.0, m_y);
266 if (m_dky) {
267 N_VDestroy_Serial(m_dky); // free derivative vector if already allocated
268 }
269 m_dky = newNVector(m_neq, m_sundials_ctx); // allocate derivative vector
270 N_VConst(0.0, m_dky);
271 // check abs tolerance array size
272 if (m_itol == CV_SV && m_nabs < m_neq) {
273 throw CanteraError("CVodesIntegrator::initialize",
274 "not enough absolute tolerance values specified.");
275 }
276
277 func.getState(NV_DATA_S(m_y));
278
279 if (m_cvode_mem) {
280 CVodeFree(&m_cvode_mem);
281 }
282
283 //! Specify the method and the iteration type. Cantera Defaults:
284 //! CV_BDF - Use BDF methods
285 //! CV_NEWTON - use Newton's method
286 #if SUNDIALS_VERSION_MAJOR < 6
287 m_cvode_mem = CVodeCreate(m_method);
288 #else
289 m_cvode_mem = CVodeCreate(m_method, m_sundials_ctx.get());
290 #endif
291 if (!m_cvode_mem) {
292 throw CanteraError("CVodesIntegrator::initialize",
293 "CVodeCreate failed.");
294 }
295
296 int flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, m_y);
297 if (flag != CV_SUCCESS) {
298 if (flag == CV_MEM_FAIL) {
299 throw CanteraError("CVodesIntegrator::initialize",
300 "Memory allocation failed.");
301 } else if (flag == CV_ILL_INPUT) {
302 throw CanteraError("CVodesIntegrator::initialize",
303 "Illegal value for CVodeInit input argument.");
304 } else {
305 throw CanteraError("CVodesIntegrator::initialize",
306 "CVodeInit failed.");
307 }
308 }
309 #if SUNDIALS_VERSION_MAJOR >= 7
310 flag = SUNContext_PushErrHandler(m_sundials_ctx.get(), &sundials_err, this);
311 #else
312 flag = CVodeSetErrHandlerFn(m_cvode_mem, &cvodes_err, this);
313 #endif
314
315 if (m_itol == CV_SV) {
316 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
317 checkError(flag, "initialize", "CVodeSVtolerances");
318 } else {
319 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
320 checkError(flag, "initialize", "CVodeSStolerances");
321 }
322
323 flag = CVodeSetUserData(m_cvode_mem, &func);
324 checkError(flag, "initialize", "CVodeSetUserData");
325
326 if (func.nparams() > 0) {
327 sensInit(t0, func);
328 flag = CVodeSetSensParams(m_cvode_mem, func.m_sens_params.data(),
329 func.m_paramScales.data(), NULL);
330 checkError(flag, "initialize", "CVodeSetSensParams");
331 }
332 applyOptions();
333}
334
335void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
336{
337 m_t0 = t0;
338 m_time = t0;
339 m_tInteg = t0;
340 func.getState(NV_DATA_S(m_y));
341 m_func = &func;
342 func.clearErrors();
343 // reinitialize preconditioner if applied
344 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
345 m_preconditioner->initialize(m_neq);
346 }
347 int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
348 checkError(result, "reinitialize", "CVodeReInit");
349 applyOptions();
350}
351
353{
354 if (m_type == "DENSE") {
355 sd_size_t N = static_cast<sd_size_t>(m_neq);
356 SUNLinSolFree((SUNLinearSolver) m_linsol);
357 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
358 #if SUNDIALS_VERSION_MAJOR >= 6
359 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
360 #else
361 m_linsol_matrix = SUNDenseMatrix(N, N);
362 #endif
363 if (m_linsol_matrix == nullptr) {
364 throw CanteraError("CVodesIntegrator::applyOptions",
365 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
366 }
367 int flag;
368 #if SUNDIALS_VERSION_MAJOR >= 6
369 #if CT_SUNDIALS_USE_LAPACK
370 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
371 m_sundials_ctx.get());
372 #else
373 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
374 m_sundials_ctx.get());
375 #endif
376 flag = CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
377 (SUNMatrix) m_linsol_matrix);
378 #else
379 #if CT_SUNDIALS_USE_LAPACK
380 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
381 #else
382 m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
383 #endif
384 flag = CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
385 (SUNMatrix) m_linsol_matrix);
386 #endif
387 if (m_linsol == nullptr) {
388 throw CanteraError("CVodesIntegrator::applyOptions",
389 "Error creating Sundials dense linear solver object");
390 } else if (flag != CV_SUCCESS) {
391 throw CanteraError("CVodesIntegrator::applyOptions",
392 "Error connecting linear solver to CVODES. "
393 "Sundials error code: {}", flag);
394 }
395
396 // throw preconditioner error for DENSE + NOJAC
397 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
398 throw CanteraError("CVodesIntegrator::applyOptions",
399 "Preconditioning is not available with the specified problem type.");
400 }
401 } else if (m_type == "DIAG") {
402 CVDiag(m_cvode_mem);
403 // throw preconditioner error for DIAG
404 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
405 throw CanteraError("CVodesIntegrator::applyOptions",
406 "Preconditioning is not available with the specified problem type.");
407 }
408 } else if (m_type == "GMRES") {
409 SUNLinSolFree((SUNLinearSolver) m_linsol);
410 #if SUNDIALS_VERSION_MAJOR >= 6
411 m_linsol = SUNLinSol_SPGMR(m_y, SUN_PREC_NONE, 0, m_sundials_ctx.get());
412 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol, nullptr);
413 #else
414 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
415 CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
416 #endif
417 // set preconditioner if used
418 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
419 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver) m_linsol,
420 static_cast<int>(m_prec_side));
421 CVodeSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
422 cvodes_prec_solve);
423 }
424 } else if (m_type == "BAND") {
425 sd_size_t N = static_cast<sd_size_t>(m_neq);
426 sd_size_t nu = m_mupper;
427 sd_size_t nl = m_mlower;
428 SUNLinSolFree((SUNLinearSolver) m_linsol);
429 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
430 #if SUNDIALS_VERSION_MAJOR >= 6
431 m_linsol_matrix = SUNBandMatrix(N, nu, nl, m_sundials_ctx.get());
432 #else
433 m_linsol_matrix = SUNBandMatrix(N, nu, nl);
434 #endif
435 if (m_linsol_matrix == nullptr) {
436 throw CanteraError("CVodesIntegrator::applyOptions",
437 "Unable to create SUNBandMatrix of size {} with bandwidths "
438 "{} and {}", N, nu, nl);
439 }
440 #if SUNDIALS_VERSION_MAJOR >= 6
441 #if CT_SUNDIALS_USE_LAPACK
442 m_linsol = SUNLinSol_LapackBand(m_y, (SUNMatrix) m_linsol_matrix,
443 m_sundials_ctx.get());
444 #else
445 m_linsol = SUNLinSol_Band(m_y, (SUNMatrix) m_linsol_matrix,
446 m_sundials_ctx.get());
447 #endif
448 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
449 (SUNMatrix) m_linsol_matrix);
450 #else
451 #if CT_SUNDIALS_USE_LAPACK
452 m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
453 #else
454 m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
455 #endif
456 CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
457 (SUNMatrix) m_linsol_matrix);
458 #endif
459 } else {
460 throw CanteraError("CVodesIntegrator::applyOptions",
461 "unsupported linear solver flag '{}'", m_type);
462 }
463
464 if (m_maxord > 0) {
465 int flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
466 checkError(flag, "applyOptions", "CVodeSetMaxOrd");
467 }
468 if (m_maxsteps > 0) {
469 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
470 }
471 if (m_hmax > 0) {
472 CVodeSetMaxStep(m_cvode_mem, m_hmax);
473 }
474 if (m_hmin > 0) {
475 CVodeSetMinStep(m_cvode_mem, m_hmin);
476 }
477 if (m_maxErrTestFails > 0) {
478 CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
479 }
480}
481
483{
484 if (tout == m_time) {
485 return;
486 } else if (tout < m_time) {
487 throw CanteraError("CVodesIntegrator::integrate",
488 "Cannot integrate backwards in time.\n"
489 "Requested time {} < current time {}",
490 tout, m_time);
491 }
492 int nsteps = 0;
493 while (m_tInteg < tout) {
494 if (nsteps >= m_maxsteps) {
495 string f_errs = m_func->getErrors();
496 if (!f_errs.empty()) {
497 f_errs = "\nExceptions caught during RHS evaluation:\n" + f_errs;
498 }
499 throw CanteraError("CVodesIntegrator::integrate",
500 "Maximum number of timesteps ({}) taken without reaching output "
501 "time ({}).\nCurrent integrator time: {}{}",
502 nsteps, tout, m_tInteg, f_errs);
503 }
504 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
505 if (flag != CV_SUCCESS) {
506 string f_errs = m_func->getErrors();
507 if (!f_errs.empty()) {
508 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
509 }
510 throw CanteraError("CVodesIntegrator::integrate",
511 "CVodes error encountered. Error code: {}\n{}\n"
512 "{}"
513 "Components with largest weighted error estimates:\n{}",
514 flag, m_error_message, f_errs, getErrorInfo(10));
515 }
516 nsteps++;
517 }
518 int flag = CVodeGetDky(m_cvode_mem, tout, 0, m_y);
519 checkError(flag, "integrate", "CVodeGetDky");
520 m_time = tout;
521 m_sens_ok = false;
522}
523
524double CVodesIntegrator::step(double tout)
525{
526 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
527 if (flag != CV_SUCCESS) {
528 string f_errs = m_func->getErrors();
529 if (!f_errs.empty()) {
530 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
531 }
532 throw CanteraError("CVodesIntegrator::step",
533 "CVodes error encountered. Error code: {}\n{}\n"
534 "{}"
535 "Components with largest weighted error estimates:\n{}",
536 flag, f_errs, m_error_message, getErrorInfo(10));
537
538 }
539 m_sens_ok = false;
541 return m_time;
542}
543
544double* CVodesIntegrator::derivative(double tout, int n)
545{
546 int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
547 checkError(flag, "derivative", "CVodeGetDky");
548 return NV_DATA_S(m_dky);
549}
550
552{
553 int ord;
554 CVodeGetLastOrder(m_cvode_mem, &ord);
555 return ord;
556}
557
559{
560 long int ne;
561 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
562 return ne;
563}
564
566{
567 // AnyMap to return stats
568 AnyMap stats;
569
570 // long int linear solver stats provided by CVodes
571 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
572 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
573 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
574 nonlinConvFails = 0, orderReductions = 0;
575 int lastOrder = 0;
576;
577
578 CVodeGetNumSteps(m_cvode_mem, &steps);
579 CVodeGetNumRhsEvals(m_cvode_mem, &rhsEvals);
580 CVodeGetNonlinSolvStats(m_cvode_mem, &nonlinIters, &nonlinConvFails);
581 CVodeGetNumErrTestFails(m_cvode_mem, &errTestFails);
582 CVodeGetLastOrder(m_cvode_mem, &lastOrder);
583 CVodeGetNumStabLimOrderReds(m_cvode_mem, &orderReductions);
584 CVodeGetNumJacEvals(m_cvode_mem, &jacEvals);
585 CVodeGetNumLinRhsEvals(m_cvode_mem, &linRhsEvals);
586 CVodeGetNumLinSolvSetups(m_cvode_mem, &linSetup);
587 CVodeGetNumLinIters(m_cvode_mem, &linIters);
588 CVodeGetNumLinConvFails(m_cvode_mem, &linConvFails);
589 CVodeGetNumPrecEvals(m_cvode_mem, &precEvals);
590 CVodeGetNumPrecSolves(m_cvode_mem, &precSolves);
591 CVodeGetNumJTSetupEvals(m_cvode_mem, &jtSetupEvals);
592 CVodeGetNumJtimesEvals(m_cvode_mem, &jTimesEvals);
593
594 #if SUNDIALS_VERSION_MAJOR >= 7 || (SUNDIALS_VERSION_MAJOR == 6 && SUNDIALS_VERSION_MINOR >= 2)
595 long int stepSolveFails = 0;
596 CVodeGetNumStepSolveFails(m_cvode_mem, &stepSolveFails);
597 stats["step_solve_fails"] = stepSolveFails;
598 #endif
599
600 stats["steps"] = steps;
601 stats["rhs_evals"] = rhsEvals;
602 stats["nonlinear_iters"] = nonlinIters;
603 stats["nonlinear_conv_fails"] = nonlinConvFails;
604 stats["err_test_fails"] = errTestFails;
605 stats["last_order"] = lastOrder;
606 stats["stab_order_reductions"] = orderReductions;
607
608 stats["jac_evals"] = jacEvals;
609 stats["lin_solve_setups"] = linSetup;
610 stats["lin_rhs_evals"] = linRhsEvals;
611 stats["lin_iters"] = linIters;
612 stats["lin_conv_fails"] = linConvFails;
613 stats["prec_evals"] = precEvals;
614 stats["prec_solves"] = precSolves;
615 stats["jt_vec_setup_evals"] = jtSetupEvals;
616 stats["jt_vec_prod_evals"] = jTimesEvals;
617 return stats;
618}
619
620double CVodesIntegrator::sensitivity(size_t k, size_t p)
621{
622 if (m_time == m_t0) {
623 // calls to CVodeGetSens are only allowed after a successful time step.
624 return 0.0;
625 }
626 if (!m_sens_ok && m_np) {
627 int flag = CVodeGetSensDky(m_cvode_mem, m_time, 0, m_yS);
628 checkError(flag, "sensitivity", "CVodeGetSens");
629 m_sens_ok = true;
630 }
631
632 if (k >= m_neq) {
633 throw CanteraError("CVodesIntegrator::sensitivity",
634 "sensitivity: k out of range ({})", k);
635 }
636 if (p >= m_np) {
637 throw CanteraError("CVodesIntegrator::sensitivity",
638 "sensitivity: p out of range ({})", p);
639 }
640 return NV_Ith_S(m_yS[p],k);
641}
642
644{
645 N_Vector errs = newNVector(m_neq, m_sundials_ctx);
646 N_Vector errw = newNVector(m_neq, m_sundials_ctx);
647 CVodeGetErrWeights(m_cvode_mem, errw);
648 CVodeGetEstLocalErrors(m_cvode_mem, errs);
649
650 vector<tuple<double, double, size_t>> weightedErrors;
651 for (size_t i=0; i<m_neq; i++) {
652 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
653 weightedErrors.emplace_back(-abs(err), err, i);
654 }
655 N_VDestroy(errs);
656 N_VDestroy(errw);
657
658 N = std::min(N, static_cast<int>(m_neq));
659 sort(weightedErrors.begin(), weightedErrors.end());
660 fmt::memory_buffer s;
661 for (int i=0; i<N; i++) {
662 fmt_append(s, "{}: {}\n",
663 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
664 }
665 return to_string(s);
666}
667
668void CVodesIntegrator::checkError(long flag, const string& ctMethod,
669 const string& cvodesMethod) const
670{
671 if (flag == CV_SUCCESS) {
672 return;
673 } else if (flag == CV_MEM_NULL) {
674 throw CanteraError("CVodesIntegrator::" + ctMethod,
675 "CVODES integrator is not initialized");
676 } else {
677 const char* flagname = CVodeGetReturnFlagName(flag);
678 throw CanteraError("CVodesIntegrator::" + ctMethod,
679 "{} returned error code {} ({}):\n{}",
680 cvodesMethod, flag, flagname, m_error_message);
681 }
682}
683
684}
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:432
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 * 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
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
Definition FuncEval.h:184
int preconditioner_setup_nothrow(double t, double *y, double gamma)
Preconditioner setup that doesn't throw an error but returns a CVODES flag.
Definition FuncEval.cpp:80
void clearErrors()
Clear any previously-stored suppressed errors.
Definition FuncEval.h:174
int evalNoThrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
Definition FuncEval.cpp:7
virtual void updatePreconditioner(double gamma)
Update the preconditioner based on already computed jacobian values.
Definition FuncEval.h:113
virtual size_t neq() const =0
Number of equations.
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
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
Definition FuncEval.h:142
shared_ptr< SystemJacobian > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
Definition Integrator.h:299
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
Definition Integrator.h:301
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
static int cvodes_prec_setup(sunrealtype t, N_Vector y, N_Vector ydot, sunbooleantype jok, sunbooleantype *jcurPtr, sunrealtype gamma, void *f_data)
Function called by CVodes 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
@ Adams_Method
Adams.
Definition Integrator.h:25
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.