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