Cantera  3.1.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 CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION >= 70
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 CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION >= 60
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 CT_SUNDIALS_VERSION < 40
287 m_cvode_mem = CVodeCreate(m_method, CV_NEWTON);
288 #elif CT_SUNDIALS_VERSION < 60
289 m_cvode_mem = CVodeCreate(m_method);
290 #else
291 m_cvode_mem = CVodeCreate(m_method, m_sundials_ctx.get());
292 #endif
293 if (!m_cvode_mem) {
294 throw CanteraError("CVodesIntegrator::initialize",
295 "CVodeCreate failed.");
296 }
297
298 int flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, m_y);
299 if (flag != CV_SUCCESS) {
300 if (flag == CV_MEM_FAIL) {
301 throw CanteraError("CVodesIntegrator::initialize",
302 "Memory allocation failed.");
303 } else if (flag == CV_ILL_INPUT) {
304 throw CanteraError("CVodesIntegrator::initialize",
305 "Illegal value for CVodeInit input argument.");
306 } else {
307 throw CanteraError("CVodesIntegrator::initialize",
308 "CVodeInit failed.");
309 }
310 }
311 #if CT_SUNDIALS_VERSION >= 70
312 flag = SUNContext_PushErrHandler(m_sundials_ctx.get(), &sundials_err, this);
313 #else
314 flag = CVodeSetErrHandlerFn(m_cvode_mem, &cvodes_err, this);
315 #endif
316
317 if (m_itol == CV_SV) {
318 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
319 checkError(flag, "initialize", "CVodeSVtolerances");
320 } else {
321 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
322 checkError(flag, "initialize", "CVodeSStolerances");
323 }
324
325 flag = CVodeSetUserData(m_cvode_mem, &func);
326 checkError(flag, "initialize", "CVodeSetUserData");
327
328 if (func.nparams() > 0) {
329 sensInit(t0, func);
330 flag = CVodeSetSensParams(m_cvode_mem, func.m_sens_params.data(),
331 func.m_paramScales.data(), NULL);
332 checkError(flag, "initialize", "CVodeSetSensParams");
333 }
334 applyOptions();
335}
336
337void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
338{
339 m_t0 = t0;
340 m_time = t0;
341 m_tInteg = t0;
342 func.getState(NV_DATA_S(m_y));
343 m_func = &func;
344 func.clearErrors();
345 // reinitialize preconditioner if applied
346 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
347 m_preconditioner->initialize(m_neq);
348 }
349 int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
350 checkError(result, "reinitialize", "CVodeReInit");
351 applyOptions();
352}
353
355{
356 if (m_type == "DENSE") {
357 sd_size_t N = static_cast<sd_size_t>(m_neq);
358 SUNLinSolFree((SUNLinearSolver) m_linsol);
359 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
360 #if CT_SUNDIALS_VERSION >= 60
361 m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
362 #else
363 m_linsol_matrix = SUNDenseMatrix(N, N);
364 #endif
365 if (m_linsol_matrix == nullptr) {
366 throw CanteraError("CVodesIntegrator::applyOptions",
367 "Unable to create SUNDenseMatrix of size {0} x {0}", N);
368 }
369 int flag;
370 #if CT_SUNDIALS_VERSION >= 60
371 #if CT_SUNDIALS_USE_LAPACK
372 m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
373 m_sundials_ctx.get());
374 #else
375 m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
376 m_sundials_ctx.get());
377 #endif
378 flag = CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
379 (SUNMatrix) m_linsol_matrix);
380 #else
381 #if CT_SUNDIALS_USE_LAPACK
382 m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
383 #else
384 m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
385 #endif
386 flag = CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
387 (SUNMatrix) m_linsol_matrix);
388 #endif
389 if (m_linsol == nullptr) {
390 throw CanteraError("CVodesIntegrator::applyOptions",
391 "Error creating Sundials dense linear solver object");
392 } else if (flag != CV_SUCCESS) {
393 throw CanteraError("CVodesIntegrator::applyOptions",
394 "Error connecting linear solver to CVODES. "
395 "Sundials error code: {}", flag);
396 }
397
398 // throw preconditioner error for DENSE + NOJAC
399 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
400 throw CanteraError("CVodesIntegrator::applyOptions",
401 "Preconditioning is not available with the specified problem type.");
402 }
403 } else if (m_type == "DIAG") {
404 CVDiag(m_cvode_mem);
405 // throw preconditioner error for DIAG
406 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
407 throw CanteraError("CVodesIntegrator::applyOptions",
408 "Preconditioning is not available with the specified problem type.");
409 }
410 } else if (m_type == "GMRES") {
411 #if CT_SUNDIALS_VERSION >= 60
412 m_linsol = SUNLinSol_SPGMR(m_y, SUN_PREC_NONE, 0, m_sundials_ctx.get());
413 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol, nullptr);
414 #elif CT_SUNDIALS_VERSION >= 40
415 m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
416 CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
417 # else
418 m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
419 CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
420 #endif
421 // set preconditioner if used
422 #if CT_SUNDIALS_VERSION >= 40
423 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
424 SUNLinSol_SPGMRSetPrecType((SUNLinearSolver) m_linsol,
425 static_cast<int>(m_prec_side));
426 CVodeSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
427 cvodes_prec_solve);
428 }
429 #else
430 if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
431 SUNSPGMRSetPrecType((SUNLinearSolver) m_linsol,
432 static_cast<int>(m_prec_side));
433 CVSpilsSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
434 cvodes_prec_solve);
435 }
436 #endif
437 } else if (m_type == "BAND") {
438 sd_size_t N = static_cast<sd_size_t>(m_neq);
439 sd_size_t nu = m_mupper;
440 sd_size_t nl = m_mlower;
441 SUNLinSolFree((SUNLinearSolver) m_linsol);
442 SUNMatDestroy((SUNMatrix) m_linsol_matrix);
443 #if CT_SUNDIALS_VERSION >= 60
444 m_linsol_matrix = SUNBandMatrix(N, nu, nl, m_sundials_ctx.get());
445 #elif CT_SUNDIALS_VERSION >= 40
446 m_linsol_matrix = SUNBandMatrix(N, nu, nl);
447 #else
448 m_linsol_matrix = SUNBandMatrix(N, nu, nl, nu+nl);
449 #endif
450 if (m_linsol_matrix == nullptr) {
451 throw CanteraError("CVodesIntegrator::applyOptions",
452 "Unable to create SUNBandMatrix of size {} with bandwidths "
453 "{} and {}", N, nu, nl);
454 }
455 #if CT_SUNDIALS_VERSION >= 60
456 #if CT_SUNDIALS_USE_LAPACK
457 m_linsol = SUNLinSol_LapackBand(m_y, (SUNMatrix) m_linsol_matrix,
458 m_sundials_ctx.get());
459 #else
460 m_linsol = SUNLinSol_Band(m_y, (SUNMatrix) m_linsol_matrix,
461 m_sundials_ctx.get());
462 #endif
463 CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
464 (SUNMatrix) m_linsol_matrix);
465 #else
466 #if CT_SUNDIALS_USE_LAPACK
467 m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
468 #else
469 m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
470 #endif
471 CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
472 (SUNMatrix) m_linsol_matrix);
473 #endif
474 } else {
475 throw CanteraError("CVodesIntegrator::applyOptions",
476 "unsupported linear solver flag '{}'", m_type);
477 }
478
479 if (m_maxord > 0) {
480 int flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
481 checkError(flag, "applyOptions", "CVodeSetMaxOrd");
482 }
483 if (m_maxsteps > 0) {
484 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
485 }
486 if (m_hmax > 0) {
487 CVodeSetMaxStep(m_cvode_mem, m_hmax);
488 }
489 if (m_hmin > 0) {
490 CVodeSetMinStep(m_cvode_mem, m_hmin);
491 }
492 if (m_maxErrTestFails > 0) {
493 CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
494 }
495}
496
498{
499 if (tout == m_time) {
500 return;
501 } else if (tout < m_time) {
502 throw CanteraError("CVodesIntegrator::integrate",
503 "Cannot integrate backwards in time.\n"
504 "Requested time {} < current time {}",
505 tout, m_time);
506 }
507 int nsteps = 0;
508 while (m_tInteg < tout) {
509 if (nsteps >= m_maxsteps) {
510 throw CanteraError("CVodesIntegrator::integrate",
511 "Maximum number of timesteps ({}) taken without reaching output "
512 "time ({}).\nCurrent integrator time: {}",
513 nsteps, tout, m_tInteg);
514 }
515 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
516 if (flag != CV_SUCCESS) {
517 string f_errs = m_func->getErrors();
518 if (!f_errs.empty()) {
519 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
520 }
521 throw CanteraError("CVodesIntegrator::integrate",
522 "CVodes error encountered. Error code: {}\n{}\n"
523 "{}"
524 "Components with largest weighted error estimates:\n{}",
525 flag, m_error_message, f_errs, getErrorInfo(10));
526 }
527 nsteps++;
528 }
529 int flag = CVodeGetDky(m_cvode_mem, tout, 0, m_y);
530 checkError(flag, "integrate", "CVodeGetDky");
531 m_time = tout;
532 m_sens_ok = false;
533}
534
535double CVodesIntegrator::step(double tout)
536{
537 int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
538 if (flag != CV_SUCCESS) {
539 string f_errs = m_func->getErrors();
540 if (!f_errs.empty()) {
541 f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
542 }
543 throw CanteraError("CVodesIntegrator::step",
544 "CVodes error encountered. Error code: {}\n{}\n"
545 "{}"
546 "Components with largest weighted error estimates:\n{}",
547 flag, f_errs, m_error_message, getErrorInfo(10));
548
549 }
550 m_sens_ok = false;
552 return m_time;
553}
554
555double* CVodesIntegrator::derivative(double tout, int n)
556{
557 int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
558 checkError(flag, "derivative", "CVodeGetDky");
559 return NV_DATA_S(m_dky);
560}
561
563{
564 int ord;
565 CVodeGetLastOrder(m_cvode_mem, &ord);
566 return ord;
567}
568
570{
571 long int ne;
572 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
573 return ne;
574}
575
577{
578 // AnyMap to return stats
579 AnyMap stats;
580
581 // long int linear solver stats provided by CVodes
582 long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
583 linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
584 precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
585 nonlinConvFails = 0, orderReductions = 0;
586 int lastOrder = 0;
587;
588 #if CT_SUNDIALS_VERSION >= 40
589 CVodeGetNumSteps(m_cvode_mem, &steps);
590 CVodeGetNumRhsEvals(m_cvode_mem, &rhsEvals);
591 CVodeGetNonlinSolvStats(m_cvode_mem, &nonlinIters, &nonlinConvFails);
592 CVodeGetNumErrTestFails(m_cvode_mem, &errTestFails);
593 CVodeGetLastOrder(m_cvode_mem, &lastOrder);
594 CVodeGetNumStabLimOrderReds(m_cvode_mem, &orderReductions);
595 CVodeGetNumJacEvals(m_cvode_mem, &jacEvals);
596 CVodeGetNumLinRhsEvals(m_cvode_mem, &linRhsEvals);
597 CVodeGetNumLinSolvSetups(m_cvode_mem, &linSetup);
598 CVodeGetNumLinIters(m_cvode_mem, &linIters);
599 CVodeGetNumLinConvFails(m_cvode_mem, &linConvFails);
600 CVodeGetNumPrecEvals(m_cvode_mem, &precEvals);
601 CVodeGetNumPrecSolves(m_cvode_mem, &precSolves);
602 CVodeGetNumJTSetupEvals(m_cvode_mem, &jtSetupEvals);
603 CVodeGetNumJtimesEvals(m_cvode_mem, &jTimesEvals);
604 #else
605 warn_user("CVodesIntegrator::solverStats", "Function not"
606 "supported with sundials versions less than 4.");
607 #endif
608
609 #if CT_SUNDIALS_VERION >= 60
610 long int stepSolveFails = 0;
611 CVodeGetNumStepSolveFails(m_cvode_mem, &stepSolveFails);
612 stats["step_solve_fails"] = stepSolveFails;
613 #endif
614
615 stats["steps"] = steps;
616 stats["rhs_evals"] = rhsEvals;
617 stats["nonlinear_iters"] = nonlinIters;
618 stats["nonlinear_conv_fails"] = nonlinConvFails;
619 stats["err_test_fails"] = errTestFails;
620 stats["last_order"] = lastOrder;
621 stats["stab_order_reductions"] = orderReductions;
622
623 stats["jac_evals"] = jacEvals;
624 stats["lin_solve_setups"] = linSetup;
625 stats["lin_rhs_evals"] = linRhsEvals;
626 stats["lin_iters"] = linIters;
627 stats["lin_conv_fails"] = linConvFails;
628 stats["prec_evals"] = precEvals;
629 stats["prec_solves"] = precSolves;
630 stats["jt_vec_setup_evals"] = jtSetupEvals;
631 stats["jt_vec_prod_evals"] = jTimesEvals;
632 return stats;
633}
634
635double CVodesIntegrator::sensitivity(size_t k, size_t p)
636{
637 if (m_time == m_t0) {
638 // calls to CVodeGetSens are only allowed after a successful time step.
639 return 0.0;
640 }
641 if (!m_sens_ok && m_np) {
642 int flag = CVodeGetSensDky(m_cvode_mem, m_time, 0, m_yS);
643 checkError(flag, "sensitivity", "CVodeGetSens");
644 m_sens_ok = true;
645 }
646
647 if (k >= m_neq) {
648 throw CanteraError("CVodesIntegrator::sensitivity",
649 "sensitivity: k out of range ({})", k);
650 }
651 if (p >= m_np) {
652 throw CanteraError("CVodesIntegrator::sensitivity",
653 "sensitivity: p out of range ({})", p);
654 }
655 return NV_Ith_S(m_yS[p],k);
656}
657
659{
660 N_Vector errs = newNVector(m_neq, m_sundials_ctx);
661 N_Vector errw = newNVector(m_neq, m_sundials_ctx);
662 CVodeGetErrWeights(m_cvode_mem, errw);
663 CVodeGetEstLocalErrors(m_cvode_mem, errs);
664
665 vector<tuple<double, double, size_t>> weightedErrors;
666 for (size_t i=0; i<m_neq; i++) {
667 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
668 weightedErrors.emplace_back(-abs(err), err, i);
669 }
670 N_VDestroy(errs);
671 N_VDestroy(errw);
672
673 N = std::min(N, static_cast<int>(m_neq));
674 sort(weightedErrors.begin(), weightedErrors.end());
675 fmt::memory_buffer s;
676 for (int i=0; i<N; i++) {
677 fmt_append(s, "{}: {}\n",
678 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
679 }
680 return to_string(s);
681}
682
683void CVodesIntegrator::checkError(long flag, const string& ctMethod,
684 const string& cvodesMethod) const
685{
686 if (flag == CV_SUCCESS) {
687 return;
688 } else if (flag == CV_MEM_NULL) {
689 throw CanteraError("CVodesIntegrator::" + ctMethod,
690 "CVODES integrator is not initialized");
691 } else {
692 const char* flagname = CVodeGetReturnFlagName(flag);
693 throw CanteraError("CVodesIntegrator::" + ctMethod,
694 "{} returned error code {} ({}):\n{}",
695 cvodesMethod, flag, flagname, m_error_message);
696 }
697}
698
699}
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
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< PreconditionerBase > 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.
void fmt_append(fmt::memory_buffer &b, 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
static int cvodes_rhs(sunrealtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:267
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
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.