PartMC
2.1.5
|
00001 /* Copyright (C) 2009-2011 Matthew West 00002 * Licensed under the GNU General Public License version 2 or (at your 00003 * option) any later version. See the file COPYING for details. 00004 */ 00005 00006 /** \file 00007 * \brief Interface to SUNDIALS ODE solver library for condensation. 00008 */ 00009 00010 #include <stdlib.h> 00011 #include <stdio.h> 00012 #include <cvode/cvode.h> 00013 #include <nvector/nvector_serial.h> 00014 #include <sundials/sundials_types.h> 00015 #include <cvode/cvode_impl.h> 00016 00017 /** \brief Result code indicating successful completion. 00018 */ 00019 #define PMC_CONDENSE_SOLVER_SUCCESS 0 00020 /** \brief Result code indicating failure to allocate \c y vector. 00021 */ 00022 #define PMC_CONDENSE_SOLVER_INIT_Y 1 00023 /** \brief Result code indicating failure to allocate \c abstol vector. 00024 */ 00025 #define PMC_CONDENSE_SOLVER_INIT_ABSTOL 2 00026 /** \brief Result code indicating failure to create the solver. 00027 */ 00028 #define PMC_CONDENSE_SOLVER_INIT_CVODE_MEM 3 00029 /** \brief Result code indicating failure to initialize the solver. 00030 */ 00031 #define PMC_CONDENSE_SOLVER_INIT_CVODE 4 00032 /** \brief Result code indicating failure to set tolerances. 00033 */ 00034 #define PMC_CONDENSE_SOLVER_SVTOL 5 00035 /** \brief Result code indicating failure to set maximum steps. 00036 */ 00037 #define PMC_CONDENSE_SOLVER_SET_MAX_STEPS 6 00038 /** \brief Result code indicating failure of the solver. 00039 */ 00040 #define PMC_CONDENSE_SOLVER_FAIL 7 00041 00042 static int condense_vf(realtype t, N_Vector y, N_Vector ydot, void *user_data); 00043 00044 static int condense_check_flag(void *flagvalue, char *funcname, int opt); 00045 00046 /*******************************************************/ 00047 // solver block 00048 static int condense_solver_Init(CVodeMem cv_mem); 00049 00050 static int condense_solver_Setup(CVodeMem cv_mem, int convfail, N_Vector ypred, 00051 N_Vector fpred, booleantype *jcurPtr, 00052 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3); 00053 00054 static int condense_solver_Solve(CVodeMem cv_mem, N_Vector b, N_Vector weight, 00055 N_Vector ycur, N_Vector fcur); 00056 00057 static void condense_solver_Free(CVodeMem cv_mem); 00058 /*******************************************************/ 00059 00060 /** \brief Call the ODE solver. 00061 * 00062 * \param neq The number of equations. 00063 * \param x_f A pointer to a vector of \c neq variables, giving the 00064 * initial state vector on entry and the final state vector on exit. 00065 * \param abstol_f A pointer to a vector of \c neq variables, giving 00066 * the absolute error tolerance for the corresponding state vector 00067 * component. 00068 * \param reltol_f The scalar relative tolerance. 00069 * \param t_initial_f The initial time (s). 00070 * \param t_final_f The final time (s). 00071 * \return A result code (0 is success). 00072 */ 00073 int condense_solver(int neq, double *x_f, double *abstol_f, double reltol_f, 00074 double t_initial_f, double t_final_f) 00075 { 00076 realtype reltol, t_initial, t_final, t, tout; 00077 N_Vector y, abstol; 00078 void *cvode_mem; 00079 CVodeMem cv_mem; 00080 int flag, i, pretype, maxl; 00081 realtype *y_data, *abstol_data; 00082 00083 y = abstol = NULL; 00084 cvode_mem = NULL; 00085 00086 y = N_VNew_Serial(neq); 00087 if (condense_check_flag((void *)y, "N_VNew_Serial", 0)) 00088 return PMC_CONDENSE_SOLVER_INIT_Y; 00089 00090 abstol = N_VNew_Serial(neq); 00091 if (condense_check_flag((void *)abstol, "N_VNew_Serial", 0)) 00092 return PMC_CONDENSE_SOLVER_INIT_ABSTOL; 00093 00094 y_data = NV_DATA_S(y); 00095 abstol_data = NV_DATA_S(abstol); 00096 for (i = 0; i < neq; i++) { 00097 y_data[i] = x_f[i]; 00098 abstol_data[i] = abstol_f[i]; 00099 } 00100 00101 reltol = reltol_f; 00102 t_initial = t_initial_f; 00103 t_final = t_final_f; 00104 00105 cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON); 00106 if (condense_check_flag((void *)cvode_mem, "CVodeCreate", 0)) 00107 return PMC_CONDENSE_SOLVER_INIT_CVODE_MEM; 00108 00109 flag = CVodeInit(cvode_mem, condense_vf, t_initial, y); 00110 if (condense_check_flag(&flag, "CVodeInit", 1)) 00111 return PMC_CONDENSE_SOLVER_INIT_CVODE; 00112 00113 flag = CVodeSVtolerances(cvode_mem, reltol, abstol); 00114 if (condense_check_flag(&flag, "CVodeSVtolerances", 1)) 00115 return PMC_CONDENSE_SOLVER_SVTOL; 00116 00117 flag = CVodeSetMaxNumSteps(cvode_mem, 100000); 00118 if (condense_check_flag(&flag, "CVodeSetMaxNumSteps", 1)) 00119 return PMC_CONDENSE_SOLVER_SET_MAX_STEPS; 00120 00121 /*******************************************************/ 00122 // dense solver 00123 //flag = CVDense(cvode_mem, neq); 00124 //if (condense_check_flag(&flag, "CVDense", 1)) return(1); 00125 /*******************************************************/ 00126 00127 /*******************************************************/ 00128 // iterative solver 00129 //pretype = PREC_LEFT; 00130 //maxl = 0; 00131 //flag = CVSptfqmr(cvode_mem, pretype, maxl); 00132 //if (condense_check_flag(&flag, "CVSptfqmr", 1)) return(1); 00133 00134 //flag = CVSpilsSetJacTimesVecFn(cvode_mem, condense_jtimes); 00135 //if (condense_check_flag(&flag, "CVSpilsSetJacTimesVecFn", 1)) return(1); 00136 00137 //flag = CVSpilsSetPreconditioner(cvode_mem, NULL, condense_prec); 00138 //if (condense_check_flag(&flag, "CVSpilsSetPreconditioner", 1)) return(1); 00139 /*******************************************************/ 00140 00141 /*******************************************************/ 00142 // explicit solver 00143 cv_mem = (CVodeMem)cvode_mem; 00144 cv_mem->cv_linit = condense_solver_Init; 00145 cv_mem->cv_lsetup = condense_solver_Setup; 00146 cv_mem->cv_lsolve = condense_solver_Solve; 00147 cv_mem->cv_lfree = condense_solver_Free; 00148 /*******************************************************/ 00149 00150 t = t_initial; 00151 flag = CVode(cvode_mem, t_final, y, &t, CV_NORMAL); 00152 if (condense_check_flag(&flag, "CVode", 1)) 00153 return PMC_CONDENSE_SOLVER_FAIL; 00154 00155 for (i = 0; i < neq; i++) { 00156 x_f[i] = y_data[i]; 00157 } 00158 00159 N_VDestroy_Serial(y); 00160 N_VDestroy_Serial(abstol); 00161 CVodeFree(&cvode_mem); 00162 return PMC_CONDENSE_SOLVER_SUCCESS; 00163 } 00164 00165 /** \brief The ODE vector field to integrate. 00166 * 00167 * \param t The current time (s). 00168 * \param y The state vector. 00169 * \param ydot The rate of change of the state vector. 00170 * \param user_data A pointer to user-provided data. 00171 * \return A result code (0 is success). 00172 */ 00173 static int condense_vf(realtype t, N_Vector y, N_Vector ydot, void *user_data) 00174 { 00175 realtype *y_data, *ydot_data; 00176 int i, neq; 00177 double *y_f, *ydot_f; 00178 00179 neq = NV_LENGTH_S(y); 00180 y_data = NV_DATA_S(y); 00181 ydot_data = NV_DATA_S(ydot); 00182 00183 y_f = malloc(neq * sizeof(double)); 00184 ydot_f = malloc(neq * sizeof(double)); 00185 00186 for (i = 0; i < neq; i++) { 00187 y_f[i] = y_data[i]; 00188 } 00189 condense_vf_f(neq, t, y_f, ydot_f); 00190 for (i = 0; i < neq; i++) { 00191 ydot_data[i] = ydot_f[i]; 00192 } 00193 00194 free(y_f); 00195 free(ydot_f); 00196 return(0); 00197 } 00198 00199 /** \brief Check the return value from a SUNDIALS call. 00200 * 00201 * - <code>opt == 0</code> means SUNDIALS function allocates memory 00202 * so check if \c flagvalue is not a \c NULL pointer. 00203 00204 * - <code>opt == 1</code> means SUNDIALS function returns a flag so 00205 * check if the \c int pointed to by \c flagvalue has 00206 * <code>flag >= 0</code>. 00207 00208 * - <code>opt == 2</code> means function allocates memory so check 00209 * if \c flagvalue is not a \c NULL pointer. 00210 * 00211 * \param flagvalue A pointer to check (either for \c NULL, or as an 00212 * \c int pointer giving the flag value). 00213 * \param funcname A string giving the function name returning this 00214 * result code. 00215 * \param opt A flag indicating the type of check to perform. 00216 * \return A result code (0 is success). 00217 */ 00218 static int condense_check_flag(void *flagvalue, char *funcname, int opt) 00219 { 00220 int *errflag; 00221 00222 /* Check if SUNDIALS function returned NULL pointer - no memory allocated */ 00223 if (opt == 0 && flagvalue == NULL) { 00224 fprintf(stderr, "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n", 00225 funcname); 00226 return(1); } 00227 00228 /* Check if flag < 0 */ 00229 else if (opt == 1) { 00230 errflag = (int *) flagvalue; 00231 if (*errflag < 0) { 00232 fprintf(stderr, "\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n", 00233 funcname, *errflag); 00234 return(1); }} 00235 00236 /* Check if function returned NULL pointer - no memory allocated */ 00237 else if (opt == 2 && flagvalue == NULL) { 00238 fprintf(stderr, "\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n", 00239 funcname); 00240 return(1); } 00241 00242 return(0); 00243 } 00244 00245 /** \brief Initialization routine for the ODE linear equation solver. 00246 * 00247 * \param cv_mem The \c CVODE solver parameter structure. 00248 * \return A result code (0 is success). 00249 */ 00250 static int condense_solver_Init(CVodeMem cv_mem) 00251 { 00252 return(0); 00253 } 00254 00255 /** \brief Setup routine for the ODE linear equation solver. 00256 * 00257 * \param cv_mem The \c CVODE solver parameter structure. 00258 * \param convfail 00259 * \param ypred 00260 * \param fpred 00261 * \param jcurPtr 00262 * \param vtemp1 00263 * \param vtemp2 00264 * \param vtemp3 00265 * \return A result code (0 is success). 00266 */ 00267 static int condense_solver_Setup(CVodeMem cv_mem, int convfail, N_Vector ypred, 00268 N_Vector fpred, booleantype *jcurPtr, 00269 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) 00270 { 00271 return(0); 00272 } 00273 00274 /** \brief Linear solver routine for use by the ODE solver. 00275 * 00276 * Should solve the system \f$(I - \gamma J) x = b\f$, where \f$J\f$ 00277 * is the current vector field Jacobian, \f$\gamma\f$ is a given 00278 * scalar, and \f$b\f$ is a given vector. 00279 * 00280 * \param cv_mem The \c CVODE solver parameter structure. 00281 * \param b The right-hand-side of the linear system. 00282 * \param weight 00283 * \param ycur The current state vector. 00284 * \param fcur The current vector field vector. 00285 * \return A result code (0 is success). 00286 */ 00287 static int condense_solver_Solve(CVodeMem cv_mem, N_Vector b, N_Vector weight, 00288 N_Vector ycur, N_Vector fcur) 00289 { 00290 realtype *b_data, *ycur_data, *fcur_data; 00291 int i, neq; 00292 double *b_f, *ycur_f, *fcur_f; 00293 double t, gamma; 00294 00295 neq = NV_LENGTH_S(b); 00296 b_data = NV_DATA_S(b); 00297 ycur_data = NV_DATA_S(ycur); 00298 fcur_data = NV_DATA_S(fcur); 00299 00300 b_f = malloc(neq * sizeof(double)); 00301 ycur_f = malloc(neq * sizeof(double)); 00302 fcur_f = malloc(neq * sizeof(double)); 00303 00304 t = cv_mem->cv_tn; 00305 gamma = cv_mem->cv_gamma; 00306 00307 for (i = 0; i < neq; i++) { 00308 b_f[i] = b_data[i]; 00309 ycur_f[i] = ycur_data[i]; 00310 fcur_f[i] = fcur_data[i]; 00311 } 00312 condense_jac_solve_f(neq, t, ycur_f, fcur_f, b_f, gamma); 00313 for (i = 0; i < neq; i++) { 00314 b_data[i] = b_f[i]; 00315 } 00316 00317 free(b_f); 00318 free(ycur_f); 00319 free(fcur_f); 00320 return(0); 00321 } 00322 00323 /** \brief Finalization routine for the ODE linear equation solver. 00324 * 00325 * \param cv_mem The \c CVODE solver parameter structure. 00326 */ 00327 static void condense_solver_Free(CVodeMem cv_mem) 00328 { 00329 }