11 #include "../../ext/cvode/include/llnltyps.h"
12 #include "../../ext/cvode/include/llnlmath.h"
13 #include "../../ext/cvode/include/cvode.h"
14 #include "../../ext/cvode/include/cvdense.h"
15 #include "../../ext/cvode/include/cvdiag.h"
16 #include "../../ext/cvode/include/cvspgmr.h"
17 #include "../../ext/cvode/include/cvode.h"
32 static void cvode_rhs(integer N, real t, N_Vector y, N_Vector ydot,
35 double* ydata = N_VDATA(y);
36 double* ydotdata = N_VDATA(ydot);
38 f->
eval(t, ydata, ydotdata, NULL);
46 static void cvode_jac(integer N, DenseMat J, RhsFn f,
void* f_data,
47 real t, N_Vector y, N_Vector fy, N_Vector ewt, real h, real uround,
48 void* jac_data,
long int* nfePtr, N_Vector vtemp1, N_Vector vtemp2,
52 double* ydata = N_VDATA(y);
53 double* fydata = N_VDATA(fy);
54 double* ewtdata = N_VDATA(ewt);
55 double* ydot = N_VDATA(vtemp1);
62 for (j=0; j < N; j++) {
66 ydata[j] = ysave + dy;
67 dy = ydata[j] - ysave;
68 func->
eval(t, ydata, ydot, NULL);
69 for (i=0; i < N; i++) {
70 col_j[i] = (ydot[i] - fydata[i])/dy;
79 CVodeInt::CVodeInt() : m_neq(0),
95 m_ropt.resize(OPT_SIZE,0.0);
96 m_iopt =
new long[OPT_SIZE];
97 fill(m_iopt, m_iopt+OPT_SIZE,0);
100 CVodeInt::~CVodeInt()
103 CVodeFree(m_cvode_mem);
116 return N_VIth(m_y,
int(k));
127 if (m_nabs != m_neq) {
131 m_abstol = N_VNew(m_nabs, 0);
133 for (
int i=0; i<m_nabs; i++) {
134 N_VIth(m_abstol, i) = abstol[i];
174 void CVodeInt::setMaxSteps(
int nmax)
177 m_iopt[MXSTEP] = m_maxsteps;
193 m_neq = int(func.
neq());
199 m_y = N_VNew(m_neq, 0);
201 if (m_itol == 1 && m_nabs < m_neq) {
202 throw CVodeErr(
"not enough absolute tolerance values specified.");
207 m_iopt[MXSTEP] = m_maxsteps;
208 m_iopt[MAXORD] = m_maxord;
209 m_ropt[HMAX] = m_hmax;
212 CVodeFree(m_cvode_mem);
216 m_data = (
void*)&func;
219 m_cvode_mem = CVodeMalloc(m_neq,
cvode_rhs, m_t0, m_y, m_method,
220 m_iter, m_itol, &m_reltol,
221 m_abstol, m_data, NULL, 1, m_iopt,
224 m_cvode_mem = CVodeMalloc(m_neq,
cvode_rhs, m_t0, m_y, m_method,
225 m_iter, m_itol, &m_reltol,
226 &m_abstols, m_data, NULL, 1, m_iopt,
231 throw CVodeErr(
"CVodeMalloc failed.");
234 if (m_type == DENSE + NOJAC) {
235 CVDense(m_cvode_mem, NULL, NULL);
236 }
else if (m_type == DENSE + JAC) {
238 }
else if (m_type == DIAG) {
240 }
else if (m_type == GMRES) {
241 CVSpgmr(m_cvode_mem, NONE, MODIFIED_GS, 0, 0.0,
244 throw CVodeErr(
"unsupported option");
248 void CVodeInt::reinitialize(
double t0,
FuncEval& func)
254 m_iopt[MXSTEP] = m_maxsteps;
255 m_iopt[MAXORD] = m_maxord;
256 m_ropt[HMAX] = m_hmax;
259 m_data = (
void*)&func;
262 result = CVReInit(m_cvode_mem,
cvode_rhs, m_t0, m_y, m_method,
263 m_iter, m_itol, &m_reltol,
264 m_abstol, m_data, NULL, 1, m_iopt,
267 result = CVReInit(m_cvode_mem,
cvode_rhs, m_t0, m_y, m_method,
268 m_iter, m_itol, &m_reltol,
269 &m_abstols, m_data, NULL, 1, m_iopt,
274 throw CVodeErr(
"CVReInit failed.");
277 if (m_type == DENSE + NOJAC) {
278 CVDense(m_cvode_mem, NULL, NULL);
279 }
else if (m_type == DENSE + JAC) {
281 }
else if (m_type == DIAG) {
283 }
else if (m_type == GMRES) {
284 CVSpgmr(m_cvode_mem, NONE, MODIFIED_GS, 0, 0.0,
287 throw CVodeErr(
"unsupported option");
295 flag = CVode(m_cvode_mem, tout, m_y, &t, NORMAL);
296 if (flag != SUCCESS) {
297 throw CVodeErr(
" CVode error encountered. Error code: " +
int2str(flag));
305 flag = CVode(m_cvode_mem, tout, m_y, &t, ONE_STEP);
306 if (flag != SUCCESS) {
307 throw CVodeErr(
" CVode error encountered. Error code: " +
int2str(flag));
Backward Differentiation.
Exception thrown when a CVODE error is encountered.
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
static void cvode_jac(integer N, DenseMat J, RhsFn f, void *f_data, real t, N_Vector y, N_Vector fy, N_Vector ewt, real h, real uround, void *jac_data, long int *nfePtr, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
Function called by CVODE to evaluate the Jacobian matrix.
virtual void eval(double t, double *y, double *ydot, double *p)=0
Evaluate the right-hand-side function.
virtual double * solution()
The current value of the solution of the system of equations.
virtual void setMethod(MethodType t)
Set the solution method.
virtual void integrate(double tout)
Integrate the system of equations.
virtual size_t neq()=0
Number of equations.
virtual void setIterator(IterType t)
Set the linear iterator.
virtual int nEvals() const
The number of function evaluations.
virtual void setProblemType(int probtype)
Set the problem type.
static void cvode_rhs(integer N, real t, N_Vector y, N_Vector ydot, void *f_data)
Function called by CVODE to evaluate ydot given y.
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set or reset the number of equations.
virtual void initialize(double t0, FuncEval &func)
Initialize the integrator for a new problem.
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
IterType
Specifies the method used for iteration.
virtual void setMinStepSize(double hmin)
Set the minimum step size.
Contains declarations for string manipulation functions within Cantera.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Virtual base class for ODE right-hand-side function evaluators.
virtual doublereal step(double tout)
Integrate the system of equations.
MethodType
Specifies the method used to integrate the system of equations.
virtual void getInitialConditions(double t0, size_t leny, double *y)=0
Fill the solution vector with the initial conditions at initial time t0.