12 #include <cvode/cvode.h>
13 #include <nvector/nvector_serial.h>
14 #include <sundials/sundials_types.h>
15 #include <cvode/cvode_impl.h>
19 #define PMC_CONDENSE_SOLVER_SUCCESS 0
22 #define PMC_CONDENSE_SOLVER_INIT_Y 1
25 #define PMC_CONDENSE_SOLVER_INIT_ABSTOL 2
28 #define PMC_CONDENSE_SOLVER_INIT_CVODE_MEM 3
31 #define PMC_CONDENSE_SOLVER_INIT_CVODE 4
34 #define PMC_CONDENSE_SOLVER_SVTOL 5
37 #define PMC_CONDENSE_SOLVER_SET_MAX_STEPS 6
40 #define PMC_CONDENSE_SOLVER_FAIL 7
42 static int condense_vf(realtype t, N_Vector y, N_Vector ydot,
void *user_data);
51 N_Vector fpred, booleantype *jcurPtr,
52 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3);
55 N_Vector ycur, N_Vector fcur);
74 double t_initial_f,
double t_final_f)
76 realtype reltol, t_initial, t_final, t, tout;
80 int flag, i, pretype, maxl;
81 realtype *y_data, *abstol_data;
86 y = N_VNew_Serial(neq);
90 abstol = N_VNew_Serial(neq);
94 y_data = NV_DATA_S(y);
95 abstol_data = NV_DATA_S(abstol);
96 for (i = 0; i < neq; i++) {
98 abstol_data[i] = abstol_f[i];
102 t_initial = t_initial_f;
105 cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON);
109 flag = CVodeInit(cvode_mem,
condense_vf, t_initial, y);
113 flag = CVodeSVtolerances(cvode_mem, reltol, abstol);
117 flag = CVodeSetMaxNumSteps(cvode_mem, 100000);
143 cv_mem = (CVodeMem)cvode_mem;
151 flag = CVode(cvode_mem, t_final, y, &t, CV_NORMAL);
155 for (i = 0; i < neq; i++) {
159 N_VDestroy_Serial(y);
160 N_VDestroy_Serial(abstol);
161 CVodeFree(&cvode_mem);
173 static int condense_vf(realtype t, N_Vector y, N_Vector ydot,
void *user_data)
175 realtype *y_data, *ydot_data;
177 double *y_f, *ydot_f;
179 neq = NV_LENGTH_S(y);
180 y_data = NV_DATA_S(y);
181 ydot_data = NV_DATA_S(ydot);
183 y_f = malloc(neq *
sizeof(
double));
184 ydot_f = malloc(neq *
sizeof(
double));
186 for (i = 0; i < neq; i++) {
189 condense_vf_f(neq, t, y_f, ydot_f);
190 for (i = 0; i < neq; i++) {
191 ydot_data[i] = ydot_f[i];
223 if (opt == 0 && flagvalue == NULL) {
225 "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n",
231 errflag = (
int *) flagvalue;
233 fprintf(stderr,
"\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
238 else if (opt == 2 && flagvalue == NULL) {
239 fprintf(stderr,
"\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
269 N_Vector fpred, booleantype *jcurPtr,
270 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
289 N_Vector ycur, N_Vector fcur)
291 realtype *b_data, *ycur_data, *fcur_data;
293 double *b_f, *ycur_f, *fcur_f;
296 neq = NV_LENGTH_S(b);
297 b_data = NV_DATA_S(b);
298 ycur_data = NV_DATA_S(ycur);
299 fcur_data = NV_DATA_S(fcur);
301 b_f = malloc(neq *
sizeof(
double));
302 ycur_f = malloc(neq *
sizeof(
double));
303 fcur_f = malloc(neq *
sizeof(
double));
306 gamma = cv_mem->cv_gamma;
308 for (i = 0; i < neq; i++) {
310 ycur_f[i] = ycur_data[i];
311 fcur_f[i] = fcur_data[i];
313 condense_jac_solve_f(neq, t, ycur_f, fcur_f, b_f, gamma);
314 for (i = 0; i < neq; i++) {
#define PMC_CONDENSE_SOLVER_INIT_Y
Result code indicating failure to allocate y vector.
#define PMC_CONDENSE_SOLVER_SET_MAX_STEPS
Result code indicating failure to set maximum steps.
static int condense_vf(realtype t, N_Vector y, N_Vector ydot, void *user_data)
The ODE vector field to integrate.
#define PMC_CONDENSE_SOLVER_SUCCESS
Result code indicating successful completion.
#define PMC_CONDENSE_SOLVER_SVTOL
Result code indicating failure to set tolerances.
int condense_solver(int neq, double *x_f, double *abstol_f, double reltol_f, double t_initial_f, double t_final_f)
Call the ODE solver.
static int condense_solver_Setup(CVodeMem cv_mem, int convfail, N_Vector ypred, N_Vector fpred, booleantype *jcurPtr, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
Setup routine for the ODE linear equation solver.
static int condense_check_flag(void *flagvalue, char *funcname, int opt)
Check the return value from a SUNDIALS call.
#define PMC_CONDENSE_SOLVER_INIT_CVODE_MEM
Result code indicating failure to create the solver.
#define PMC_CONDENSE_SOLVER_INIT_CVODE
Result code indicating failure to initialize the solver.
static int condense_solver_Solve(CVodeMem cv_mem, N_Vector b, N_Vector weight, N_Vector ycur, N_Vector fcur)
Linear solver routine for use by the ODE solver.
#define PMC_CONDENSE_SOLVER_INIT_ABSTOL
Result code indicating failure to allocate abstol vector.
static void condense_solver_Free(CVodeMem cv_mem)
Finalization routine for the ODE linear equation solver.
static int condense_solver_Init(CVodeMem cv_mem)
Initialization routine for the ODE linear equation solver.
#define PMC_CONDENSE_SOLVER_FAIL
Result code indicating failure of the solver.