7 #include "cantera/base/config.h"
15 #if SUNDIALS_VERSION == 22
17 #include "sundials_types.h"
18 #include "sundials_math.h"
20 #include "cvodes_dense.h"
21 #include "cvodes_diag.h"
22 #include "cvodes_spgmr.h"
23 #include "cvodes_band.h"
24 #include "nvector_serial.h"
28 #if SUNDIALS_VERSION >= 23
29 #include "sundials/sundials_types.h"
30 #include "sundials/sundials_math.h"
31 #include "sundials/sundials_nvector.h"
32 #include "nvector/nvector_serial.h"
33 #include "cvodes/cvodes.h"
34 #include "cvodes/cvodes_dense.h"
35 #include "cvodes/cvodes_diag.h"
36 #include "cvodes/cvodes_spgmr.h"
37 #include "cvodes/cvodes_band.h"
40 #error unsupported Sundials version!
43 #if SUNDIALS_VERSION >= 24
50 inline static N_Vector nv(
void* x)
52 return reinterpret_cast<N_Vector
>(x);
61 FuncData(FuncEval* f,
int npar = 0) {
62 m_pars.resize(npar, 1.0);
65 virtual ~FuncData() {}
83 static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot,
86 double* ydata = NV_DATA_S(y);
87 double* ydotdata = NV_DATA_S(ydot);
88 Cantera::FuncData* d = (Cantera::FuncData*)f_data;
91 if (d->m_pars.size() == 0) {
92 f->
eval(t, ydata, ydotdata, NULL);
114 CVodesIntegrator::CVodesIntegrator() :
127 m_reltolsens(1.0e-5),
128 m_abstolsens(1.0e-4),
134 m_mupper(0), m_mlower(0)
143 CVodesIntegrator::~CVodesIntegrator()
147 CVodeSensFree(m_cvode_mem);
149 CVodeFree(&m_cvode_mem);
152 N_VDestroy_Serial(nv(m_y));
155 N_VDestroy_Serial(nv(m_abstol));
162 double& CVodesIntegrator::solution(
size_t k)
164 return NV_Ith_S(nv(m_y),k);
167 double* CVodesIntegrator::solution()
169 return NV_DATA_S(nv(m_y));
172 void CVodesIntegrator::setTolerances(
double reltol,
size_t n,
double* abstol)
178 N_VDestroy_Serial(nv(m_abstol));
180 m_abstol =
reinterpret_cast<void*
>(N_VNew_Serial(n));
182 for (
size_t i=0; i<n; i++) {
183 NV_Ith_S(nv(m_abstol), i) = abstol[i];
188 void CVodesIntegrator::setTolerances(
double reltol,
double abstol)
195 void CVodesIntegrator::setSensitivityTolerances(
double reltol,
double abstol)
197 m_reltolsens = reltol;
198 m_abstolsens = abstol;
201 void CVodesIntegrator::setProblemType(
int probtype)
206 void CVodesIntegrator::setMethod(
MethodType t)
213 throw CVodesErr(
"unknown method");
217 void CVodesIntegrator::setMaxStepSize(doublereal hmax)
221 CVodeSetMaxStep(m_cvode_mem, hmax);
225 void CVodesIntegrator::setMinStepSize(doublereal hmin)
229 CVodeSetMinStep(m_cvode_mem, hmin);
233 void CVodesIntegrator::setMaxSteps(
int nmax)
237 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
241 void CVodesIntegrator::setIterator(
IterType t)
246 m_iter = CV_FUNCTIONAL;
248 throw CVodesErr(
"unknown iterator");
252 void CVodesIntegrator::sensInit(
double t0, FuncEval& func)
254 m_np = func.nparams();
255 size_t nv = func.neq();
259 y = N_VNew_Serial(nv);
260 m_yS = N_VCloneVectorArray_Serial(m_np, y);
261 for (
size_t n = 0; n < m_np; n++) {
262 data = NV_DATA_S(m_yS[n]);
263 for (
size_t j = 0; j < nv; j++) {
270 #if SUNDIALS_VERSION <= 23
271 flag = CVodeSensMalloc(m_cvode_mem, m_np, CV_STAGGERED, m_yS);
272 if (flag != CV_SUCCESS) {
273 throw CVodesErr(
"Error in CVodeSensMalloc");
275 vector_fp atol(m_np, m_abstolsens);
276 double rtol = m_reltolsens;
277 flag = CVodeSetSensTolerances(m_cvode_mem, CV_SS, rtol,
DATA_PTR(atol));
278 #elif SUNDIALS_VERSION >= 24
279 flag = CVodeSensInit(m_cvode_mem, m_np, CV_STAGGERED,
280 CVSensRhsFn(0), m_yS);
282 if (flag != CV_SUCCESS) {
283 throw CVodesErr(
"Error in CVodeSensMalloc");
285 vector_fp atol(m_np, m_abstolsens);
286 double rtol = m_reltolsens;
287 flag = CVodeSensSStolerances(m_cvode_mem, rtol,
DATA_PTR(atol));
292 void CVodesIntegrator::initialize(
double t0, FuncEval& func)
298 N_VDestroy_Serial(nv(m_y));
300 m_y =
reinterpret_cast<void*
>(N_VNew_Serial(m_neq));
301 for (
size_t i = 0; i < m_neq; i++) {
302 NV_Ith_S(nv(m_y), i) = 0.0;
305 if (m_itol == CV_SV && m_nabs < m_neq) {
306 throw CVodesErr(
"not enough absolute tolerance values specified.");
309 func.getInitialConditions(m_t0, m_neq, NV_DATA_S(nv(m_y)));
312 CVodeFree(&m_cvode_mem);
321 m_cvode_mem = CVodeCreate(m_method, m_iter);
323 throw CVodesErr(
"CVodeCreate failed.");
327 #if SUNDIALS_VERSION <= 23
328 if (m_itol == CV_SV) {
330 flag = CVodeMalloc(m_cvode_mem,
cvodes_rhs, m_t0, nv(m_y), m_itol,
331 m_reltol, nv(m_abstol));
334 flag = CVodeMalloc(m_cvode_mem,
cvodes_rhs, m_t0, nv(m_y), m_itol,
335 m_reltol, &m_abstols);
337 if (flag != CV_SUCCESS) {
338 if (flag == CV_MEM_FAIL) {
339 throw CVodesErr(
"Memory allocation failed.");
340 }
else if (flag == CV_ILL_INPUT) {
341 throw CVodesErr(
"Illegal value for CVodeMalloc input argument.");
343 throw CVodesErr(
"CVodeMalloc failed.");
346 #elif SUNDIALS_VERSION >= 24
348 flag = CVodeInit(m_cvode_mem,
cvodes_rhs, m_t0, nv(m_y));
349 if (flag != CV_SUCCESS) {
350 if (flag == CV_MEM_FAIL) {
351 throw CVodesErr(
"Memory allocation failed.");
352 }
else if (flag == CV_ILL_INPUT) {
353 throw CVodesErr(
"Illegal value for CVodeInit input argument.");
355 throw CVodesErr(
"CVodeInit failed.");
359 if (m_itol == CV_SV) {
360 flag = CVodeSVtolerances(m_cvode_mem, m_reltol, nv(m_abstol));
362 flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
364 if (flag != CV_SUCCESS) {
365 if (flag == CV_MEM_FAIL) {
366 throw CVodesErr(
"Memory allocation failed.");
367 }
else if (flag == CV_ILL_INPUT) {
368 throw CVodesErr(
"Illegal value for CVodeInit input argument.");
370 throw CVodesErr(
"CVodeInit failed.");
375 if (m_type == DENSE + NOJAC) {
377 CVDense(m_cvode_mem, N);
378 }
else if (m_type == DIAG) {
380 }
else if (m_type == GMRES) {
381 CVSpgmr(m_cvode_mem, PREC_NONE, 0);
382 }
else if (m_type == BAND + NOJAC) {
384 long int nu = m_mupper;
385 long int nl = m_mlower;
386 CVBand(m_cvode_mem, N, nu, nl);
388 throw CVodesErr(
"unsupported option");
393 m_fdata =
new FuncData(&func, func.nparams());
396 #if SUNDIALS_VERSION <= 23
397 flag = CVodeSetFdata(m_cvode_mem, (
void*)m_fdata);
398 if (flag != CV_SUCCESS) {
399 throw CVodesErr(
"CVodeSetFdata failed.");
401 #elif SUNDIALS_VERSION >= 24
402 flag = CVodeSetUserData(m_cvode_mem, (
void*)m_fdata);
403 if (flag != CV_SUCCESS) {
404 throw CVodesErr(
"CVodeSetUserData failed.");
407 if (func.nparams() > 0) {
409 flag = CVodeSetSensParams(m_cvode_mem,
DATA_PTR(m_fdata->m_pars),
415 flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
417 if (m_maxsteps > 0) {
418 flag = CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
421 flag = CVodeSetMaxStep(m_cvode_mem, m_hmax);
426 void CVodesIntegrator::reinitialize(
double t0, FuncEval& func)
430 func.getInitialConditions(m_t0, m_neq, NV_DATA_S(nv(m_y)));
439 #if SUNDIALS_VERSION <= 23
440 if (m_itol == CV_SV) {
441 result = CVodeReInit(m_cvode_mem,
cvodes_rhs, m_t0, nv(m_y),
445 result = CVodeReInit(m_cvode_mem,
cvodes_rhs, m_t0, nv(m_y),
449 if (result != CV_SUCCESS) {
450 throw CVodesErr(
"CVodeReInit failed. result = "+
int2str(result));
452 #elif SUNDIALS_VERSION >= 24
453 result = CVodeReInit(m_cvode_mem, m_t0, nv(m_y));
454 if (result != CV_SUCCESS) {
455 throw CVodesErr(
"CVodeReInit failed. result = "+
int2str(result));
459 if (m_type == DENSE + NOJAC) {
461 CVDense(m_cvode_mem, N);
462 }
else if (m_type == DIAG) {
464 }
else if (m_type == BAND + NOJAC) {
466 long int nu = m_mupper;
467 long int nl = m_mlower;
468 CVBand(m_cvode_mem, N, nu, nl);
469 }
else if (m_type == GMRES) {
470 CVSpgmr(m_cvode_mem, PREC_NONE, 0);
472 throw CVodesErr(
"unsupported option");
478 CVodeSetMaxOrd(m_cvode_mem, m_maxord);
480 if (m_maxsteps > 0) {
481 CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
484 CVodeSetMaxStep(m_cvode_mem, m_hmax);
488 void CVodesIntegrator::integrate(
double tout)
492 flag = CVode(m_cvode_mem, tout, nv(m_y), &t, CV_NORMAL);
493 if (flag != CV_SUCCESS) {
494 throw CVodesErr(
" CVodes error encountered. Error code: " +
int2str(flag));
496 #if SUNDIALS_VERSION <= 23
498 CVodeGetSens(m_cvode_mem, tout, m_yS);
500 #elif SUNDIALS_VERSION >= 24
503 CVodeGetSens(m_cvode_mem, &tretn, m_yS);
504 if (fabs(tretn - tout) > 1.0E-5) {
505 throw CVodesErr(
"Time of Sensitivities different than time of tout");
511 double CVodesIntegrator::step(
double tout)
515 flag = CVode(m_cvode_mem, tout, nv(m_y), &t, CV_ONE_STEP);
516 if (flag != CV_SUCCESS) {
517 throw CVodesErr(
" CVodes error encountered. Error code: " +
int2str(flag));
522 int CVodesIntegrator::nEvals()
const
525 CVodeGetNumRhsEvals(m_cvode_mem, &ne);
530 double CVodesIntegrator::sensitivity(
size_t k,
size_t p)
533 throw CVodesErr(
"sensitivity: k out of range ("+
int2str(p)+
")");
536 throw CVodesErr(
"sensitivity: p out of range ("+
int2str(p)+
")");
538 return NV_Ith_S(m_yS[p],k);