12 #include "sundials/sundials_types.h" 13 #include "sundials/sundials_math.h" 14 #include "sundials/sundials_nvector.h" 15 #include "nvector/nvector_serial.h" 16 #include "cvodes/cvodes.h" 17 #if SUNDIALS_USE_LAPACK 18 #include "cvodes/cvodes_lapack.h" 20 #include "cvodes/cvodes_dense.h" 21 #include "cvodes/cvodes_band.h" 23 #include "cvodes/cvodes_diag.h" 24 #include "cvodes/cvodes_spgmr.h" 29 #if SUNDIALS_VERSION < 25 30 typedef int sd_size_t;
32 typedef long int sd_size_t;
47 static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot,
void* f_data)
50 return f->
eval_nothrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
56 static void cvodes_err(
int error_code,
const char* module,
57 const char*
function,
char* msg,
void* eh_data)
65 CVodesIntegrator::CVodesIntegrator() :
88 m_mupper(0), m_mlower(0),
93 CVodesIntegrator::~CVodesIntegrator()
97 CVodeSensFree(m_cvode_mem);
99 CVodeFree(&m_cvode_mem);
102 N_VDestroy_Serial(m_y);
105 N_VDestroy_Serial(m_abstol);
108 N_VDestroyVectorArray_Serial(m_yS, static_cast<sd_size_t>(m_np));
114 return NV_Ith_S(m_y, k);
119 return NV_DATA_S(m_y);
128 N_VDestroy_Serial(m_abstol);
130 m_abstol = N_VNew_Serial(static_cast<sd_size_t>(n));
132 for (
size_t i=0; i<n; i++) {
133 NV_Ith_S(m_abstol, i) = abstol[i];
147 m_reltolsens = reltol;
148 m_abstolsens = abstol;
163 throw CanteraError(
"CVodesIntegrator::setMethod",
"unknown method");
171 CVodeSetMaxStep(m_cvode_mem, hmax);
179 CVodeSetMinStep(m_cvode_mem, hmin);
183 void CVodesIntegrator::setMaxSteps(
int nmax)
187 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
193 m_maxErrTestFails = n;
195 CVodeSetMaxErrTestFails(m_cvode_mem, n);
204 m_iter = CV_FUNCTIONAL;
206 throw CanteraError(
"CVodesIntegrator::setIterator",
"unknown iterator");
210 void CVodesIntegrator::sensInit(
double t0,
FuncEval& func)
215 N_Vector y = N_VNew_Serial(static_cast<sd_size_t>(func.
neq()));
216 m_yS = N_VCloneVectorArray_Serial(static_cast<sd_size_t>(m_np), y);
217 for (
size_t n = 0; n < m_np; n++) {
218 N_VConst(0.0, m_yS[n]);
220 N_VDestroy_Serial(y);
222 int flag = CVodeSensInit(m_cvode_mem, static_cast<sd_size_t>(m_np),
223 CV_STAGGERED, CVSensRhsFn(0), m_yS);
225 if (flag != CV_SUCCESS) {
226 throw CanteraError(
"CVodesIntegrator::sensInit",
"Error in CVodeSensInit");
229 for (
size_t n = 0; n < m_np; n++) {
234 flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
246 N_VDestroy_Serial(m_y);
248 m_y = N_VNew_Serial(static_cast<sd_size_t>(m_neq));
251 if (m_itol == CV_SV && m_nabs < m_neq) {
253 "not enough absolute tolerance values specified.");
259 CVodeFree(&m_cvode_mem);
265 m_cvode_mem = CVodeCreate(m_method, m_iter);
268 "CVodeCreate failed.");
271 int flag = CVodeInit(m_cvode_mem,
cvodes_rhs, m_t0, m_y);
272 if (flag != CV_SUCCESS) {
273 if (flag == CV_MEM_FAIL) {
275 "Memory allocation failed.");
276 }
else if (flag == CV_ILL_INPUT) {
278 "Illegal value for CVodeInit input argument.");
281 "CVodeInit failed.");
284 CVodeSetErrHandlerFn(m_cvode_mem, &
cvodes_err,
this);
286 if (m_itol == CV_SV) {
287 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
289 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
291 if (flag != CV_SUCCESS) {
292 if (flag == CV_MEM_FAIL) {
294 "Memory allocation failed.");
295 }
else if (flag == CV_ILL_INPUT) {
297 "Illegal value for CVodeInit input argument.");
300 "CVodeInit failed.");
304 flag = CVodeSetUserData(m_cvode_mem, &func);
305 if (flag != CV_SUCCESS) {
307 "CVodeSetUserData failed.");
311 flag = CVodeSetSensParams(m_cvode_mem, func.
m_sens_params.data(),
313 if (flag != CV_SUCCESS) {
315 "CVodeSetSensParams failed.");
321 void CVodesIntegrator::reinitialize(
double t0,
FuncEval& func)
329 int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
330 if (result != CV_SUCCESS) {
332 "CVodeReInit failed. result = {}", result);
339 if (m_type == DENSE + NOJAC) {
340 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
341 #if SUNDIALS_USE_LAPACK 342 CVLapackDense(m_cvode_mem, N);
344 CVDense(m_cvode_mem, N);
346 }
else if (m_type == DIAG) {
348 }
else if (m_type == GMRES) {
349 CVSpgmr(m_cvode_mem, PREC_NONE, 0);
350 }
else if (m_type == BAND + NOJAC) {
351 sd_size_t N =
static_cast<sd_size_t
>(m_neq);
352 long int nu = m_mupper;
353 long int nl = m_mlower;
354 #if SUNDIALS_USE_LAPACK 355 CVLapackBand(m_cvode_mem, N, nu, nl);
357 CVBand(m_cvode_mem, N, nu, nl);
361 "unsupported option");
365 CVodeSetMaxOrd(m_cvode_mem, m_maxord);
367 if (m_maxsteps > 0) {
368 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
371 CVodeSetMaxStep(m_cvode_mem, m_hmax);
374 CVodeSetMinStep(m_cvode_mem, m_hmin);
376 if (m_maxErrTestFails > 0) {
377 CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
386 int flag = CVode(m_cvode_mem, tout, m_y, &
m_time, CV_NORMAL);
387 if (flag != CV_SUCCESS) {
389 if (!f_errs.empty()) {
390 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
393 "CVodes error encountered. Error code: {}\n{}\n" 395 "Components with largest weighted error estimates:\n{}",
403 int flag = CVode(m_cvode_mem, tout, m_y, &
m_time, CV_ONE_STEP);
404 if (flag != CV_SUCCESS) {
406 if (!f_errs.empty()) {
407 f_errs =
"Exceptions caught during RHS evaluation:\n" + f_errs;
410 "CVodes error encountered. Error code: {}\n{}\n" 412 "Components with largest weighted error estimates:\n{}",
423 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
427 double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
434 int flag = CVodeGetSens(m_cvode_mem, &
m_time, m_yS);
435 if (flag != CV_SUCCESS) {
436 throw CanteraError(
"CVodesIntegrator::sensitivity",
437 "CVodeGetSens failed. Error code: {}", flag);
443 throw CanteraError(
"CVodesIntegrator::sensitivity",
444 "sensitivity: k out of range ({})", k);
447 throw CanteraError(
"CVodesIntegrator::sensitivity",
448 "sensitivity: p out of range ({})", p);
450 return NV_Ith_S(m_yS[p],k);
455 N_Vector errs = N_VNew_Serial(static_cast<sd_size_t>(m_neq));
456 N_Vector errw = N_VNew_Serial(static_cast<sd_size_t>(m_neq));
457 CVodeGetErrWeights(m_cvode_mem, errw);
458 CVodeGetEstLocalErrors(m_cvode_mem, errs);
460 vector<tuple<double, double, size_t> > weightedErrors;
461 for (
size_t i=0; i<m_neq; i++) {
462 double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
463 weightedErrors.emplace_back(-abs(err), err, i);
468 N = std::min(N, static_cast<int>(m_neq));
469 sort(weightedErrors.begin(), weightedErrors.end());
471 for (
int i=0; i<N; i++) {
473 get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
Backward Differentiation.
virtual void setProblemType(int probtype)
Set the problem type.
void applyOptions()
Applies user-specified options to the underlying CVODES solver.
Wrapper class for 'cvodes' integrator from LLNL.
virtual int nEvals() const
The number of function evaluations.
vector_fp m_paramScales
Scaling factors for each sensitivity parameter.
void clearErrors()
Clear any previously-stored suppressed errors.
std::string getErrors() const
Return a string containing the text of any suppressed errors.
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
virtual void setMethod(MethodType t)
Set the solution method.
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
virtual void setMinStepSize(double hmin)
Set the minimum step size.
virtual void setIterator(IterType t)
Set the linear iterator.
virtual double * solution()
The current value of the solution of the system of equations.
virtual size_t neq()=0
Number of equations.
virtual void integrate(double tout)
Integrate the system of equations.
virtual void setSensitivityTolerances(double reltol, double abstol)
Set the sensitivity error tolerances.
virtual std::string getErrorInfo(int N)
Returns a string listing the weighted error estimates associated with each solution component...
Base class for exceptions thrown by Cantera classes.
vector_fp m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
virtual size_t nparams()
Number of sensitivity parameters.
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.
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set error tolerances.
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
int eval_nothrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS have been updated for at the current integrator ti...
static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
IterType
Specifies the method used for iteration.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
virtual doublereal step(double tout)
Integrate the system of equations.
Contains declarations for string manipulation functions within Cantera.
std::string m_error_message
Error message information provide by CVodes.
Virtual base class for ODE right-hand-side function evaluators.
Namespace for the Cantera kernel.
MethodType
Specifies the method used to integrate the system of equations.
double m_time
The current integrator time.
virtual void initialize(double t0, FuncEval &func)
Initialize the integrator for a new problem.