00001 /* 00002 * ----------------------------------------------------------------- 00003 * $Revision: 1.2 $ 00004 * $Date: 2006/11/29 00:05:06 $ 00005 * ----------------------------------------------------------------- 00006 * Programmer: Radu Serban @ LLNL 00007 * ----------------------------------------------------------------- 00008 * Copyright (c) 2006, The Regents of the University of California. 00009 * Produced at the Lawrence Livermore National Laboratory. 00010 * All rights reserved. 00011 * For details, see the LICENSE file. 00012 * ----------------------------------------------------------------- 00013 * Common header file for the direct linear solvers in CVODES. 00014 * 00015 * Part I contains type definitions and function prototypes for 00016 * using a CVDIRECT linear solver on forward problems (IVP 00017 * integration and/or FSA) 00018 * 00019 * Part II contains type definitions and function prototypes for 00020 * using a CVDIRECT linear solver on adjoint (backward) problems 00021 * ----------------------------------------------------------------- 00022 */ 00023 00024 #ifndef _CVSDIRECT_H 00025 #define _CVSDIRECT_H 00026 00027 #ifdef __cplusplus /* wrapper to enable C++ usage */ 00028 extern "C" { 00029 #endif 00030 00031 #include <sundials/sundials_direct.h> 00032 #include <sundials/sundials_nvector.h> 00033 00034 /* 00035 * ================================================================= 00036 * C V S D I R E C T C O N S T A N T S 00037 * ================================================================= 00038 */ 00039 00040 /* 00041 * ----------------------------------------------------------------- 00042 * CVSDIRECT return values 00043 * ----------------------------------------------------------------- 00044 */ 00045 00046 #define CVDIRECT_SUCCESS 0 00047 #define CVDIRECT_MEM_NULL -1 00048 #define CVDIRECT_LMEM_NULL -2 00049 #define CVDIRECT_ILL_INPUT -3 00050 #define CVDIRECT_MEM_FAIL -4 00051 00052 /* Additional last_flag values */ 00053 00054 #define CVDIRECT_JACFUNC_UNRECVR -5 00055 #define CVDIRECT_JACFUNC_RECVR -6 00056 00057 /* Return values for the adjoint module */ 00058 00059 #define CVDIRECT_ADJMEM_NULL -101 00060 #define CVDIRECT_LMEMB_NULL -102 00061 00062 /* 00063 * ================================================================= 00064 * PART I: F O R W A R D P R O B L E M S 00065 * ================================================================= 00066 */ 00067 00068 /* 00069 * ----------------------------------------------------------------- 00070 * FUNCTION TYPES 00071 * ----------------------------------------------------------------- 00072 */ 00073 00074 /* 00075 * ----------------------------------------------------------------- 00076 * Type: CVDlsDenseJacFn 00077 * ----------------------------------------------------------------- 00078 * 00079 * A dense Jacobian approximation function Jac must be of type 00080 * CVDlsDenseJacFn. Its parameters are: 00081 * 00082 * N is the problem size. 00083 * 00084 * Jac is the dense matrix (of type DlsMat) that will be loaded 00085 * by a CVDlsDenseJacFn with an approximation to the Jacobian 00086 * matrix J = (df_i/dy_j) at the point (t,y). 00087 * 00088 * t is the current value of the independent variable. 00089 * 00090 * y is the current value of the dependent variable vector, 00091 * namely the predicted value of y(t). 00092 * 00093 * fy is the vector f(t,y). 00094 * 00095 * jac_data is a pointer to user data - the same as the jac_data 00096 * parameter passed to CVDlsSetJacFn. 00097 * 00098 * tmp1, tmp2, and tmp3 are pointers to memory allocated for 00099 * vectors of length N which can be used by a CVDlsDenseJacFn 00100 * as temporary storage or work space. 00101 * 00102 * A CVDlsDenseJacFn should return 0 if successful, a positive 00103 * value if a recoverable error occurred, and a negative value if 00104 * an unrecoverable error occurred. 00105 * 00106 * ----------------------------------------------------------------- 00107 * 00108 * NOTE: The following are two efficient ways to load a dense Jac: 00109 * (1) (with macros - no explicit data structure references) 00110 * for (j=0; j < Neq; j++) { 00111 * col_j = DENSE_COL(Jac,j); 00112 * for (i=0; i < Neq; i++) { 00113 * generate J_ij = the (i,j)th Jacobian element 00114 * col_j[i] = J_ij; 00115 * } 00116 * } 00117 * (2) (without macros - explicit data structure references) 00118 * for (j=0; j < Neq; j++) { 00119 * col_j = (Jac->data)[j]; 00120 * for (i=0; i < Neq; i++) { 00121 * generate J_ij = the (i,j)th Jacobian element 00122 * col_j[i] = J_ij; 00123 * } 00124 * } 00125 * A third way, using the DENSE_ELEM(A,i,j) macro, is much less 00126 * efficient in general. It is only appropriate for use in small 00127 * problems in which efficiency of access is NOT a major concern. 00128 * 00129 * NOTE: If the user's Jacobian routine needs other quantities, 00130 * they are accessible as follows: hcur (the current stepsize) 00131 * and ewt (the error weight vector) are accessible through 00132 * CVodeGetCurrentStep and CVodeGetErrWeights, respectively 00133 * (see cvode.h). The unit roundoff is available as 00134 * UNIT_ROUNDOFF defined in sundials_types.h. 00135 * 00136 * ----------------------------------------------------------------- 00137 */ 00138 00139 00140 typedef int (*CVDlsDenseJacFn)(int N, realtype t, 00141 N_Vector y, N_Vector fy, 00142 DlsMat Jac, void *jac_data, 00143 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3); 00144 00145 /* 00146 * ----------------------------------------------------------------- 00147 * Type: CVDlsBandJacFn 00148 * ----------------------------------------------------------------- 00149 * 00150 * A band Jacobian approximation function Jac must have the 00151 * prototype given below. Its parameters are: 00152 * 00153 * N is the length of all vector arguments. 00154 * 00155 * mupper is the upper half-bandwidth of the approximate banded 00156 * Jacobian. This parameter is the same as the mupper parameter 00157 * passed by the user to the linear solver initialization function. 00158 * 00159 * mlower is the lower half-bandwidth of the approximate banded 00160 * Jacobian. This parameter is the same as the mlower parameter 00161 * passed by the user to the linear solver initialization function. 00162 * 00163 * t is the current value of the independent variable. 00164 * 00165 * y is the current value of the dependent variable vector, 00166 * namely the predicted value of y(t). 00167 * 00168 * fy is the vector f(t,y). 00169 * 00170 * Jac is the band matrix (of type DlsMat) that will be loaded 00171 * by a CVDlsBandJacFn with an approximation to the Jacobian matrix 00172 * Jac = (df_i/dy_j) at the point (t,y). 00173 * Three efficient ways to load J are: 00174 * 00175 * (1) (with macros - no explicit data structure references) 00176 * for (j=0; j < n; j++) { 00177 * col_j = BAND_COL(Jac,j); 00178 * for (i=j-mupper; i <= j+mlower; i++) { 00179 * generate J_ij = the (i,j)th Jacobian element 00180 * BAND_COL_ELEM(col_j,i,j) = J_ij; 00181 * } 00182 * } 00183 * 00184 * (2) (with BAND_COL macro, but without BAND_COL_ELEM macro) 00185 * for (j=0; j < n; j++) { 00186 * col_j = BAND_COL(Jac,j); 00187 * for (k=-mupper; k <= mlower; k++) { 00188 * generate J_ij = the (i,j)th Jacobian element, i=j+k 00189 * col_j[k] = J_ij; 00190 * } 00191 * } 00192 * 00193 * (3) (without macros - explicit data structure references) 00194 * offset = Jac->smu; 00195 * for (j=0; j < n; j++) { 00196 * col_j = ((Jac->data)[j])+offset; 00197 * for (k=-mupper; k <= mlower; k++) { 00198 * generate J_ij = the (i,j)th Jacobian element, i=j+k 00199 * col_j[k] = J_ij; 00200 * } 00201 * } 00202 * Caution: Jac->smu is generally NOT the same as mupper. 00203 * 00204 * The BAND_ELEM(A,i,j) macro is appropriate for use in small 00205 * problems in which efficiency of access is NOT a major concern. 00206 * 00207 * jac_data is a pointer to user data - the same as the jac_data 00208 * parameter passed to CVDlsSetJacFn. 00209 * 00210 * NOTE: If the user's Jacobian routine needs other quantities, 00211 * they are accessible as follows: hcur (the current stepsize) 00212 * and ewt (the error weight vector) are accessible through 00213 * CVodeGetCurrentStep and CVodeGetErrWeights, respectively 00214 * (see cvode.h). The unit roundoff is available as 00215 * UNIT_ROUNDOFF defined in sundials_types.h 00216 * 00217 * tmp1, tmp2, and tmp3 are pointers to memory allocated for 00218 * vectors of length N which can be used by a CVDlsBandJacFn 00219 * as temporary storage or work space. 00220 * 00221 * A CVDlsBandJacFn should return 0 if successful, a positive value 00222 * if a recoverable error occurred, and a negative value if an 00223 * unrecoverable error occurred. 00224 * ----------------------------------------------------------------- 00225 */ 00226 00227 typedef int (*CVDlsBandJacFn)(int N, int mupper, int mlower, 00228 realtype t, N_Vector y, N_Vector fy, 00229 DlsMat Jac, void *jac_data, 00230 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3); 00231 00232 /* 00233 * ----------------------------------------------------------------- 00234 * EXPORTED FUNCTIONS 00235 * ----------------------------------------------------------------- 00236 */ 00237 00238 /* 00239 * ----------------------------------------------------------------- 00240 * Optional inputs to the CVSDIRECT linear solver 00241 * ----------------------------------------------------------------- 00242 * 00243 * CVDlsSetJacFn specifies the Jacobian approximation routine to be 00244 * used. When using dense Jacobians, a user-supplied jac routine must 00245 * be of type CVDlsDenseJacFn. When using banded Jacobians, a 00246 * user-supplied jac routine must be of type CVDlsBandJacFn. 00247 * By default, a difference quotient approximation, supplied with this 00248 * solver is used. 00249 * CVDlsSetJacFn also specifies a pointer to user data which is 00250 * passed to the user's jac routine every time it is called. 00251 * 00252 * The return value of CVDlsSetJacFn is one of: 00253 * CVDIRECT_SUCCESS if successful 00254 * CVDIRECT_MEM_NULL if the CVODES memory was NULL 00255 * CVDIRECT_LMEM_NULL if the linear solver memory was NULL 00256 * ----------------------------------------------------------------- 00257 */ 00258 00259 SUNDIALS_EXPORT int CVDlsSetJacFn(void *cvode_mem, void *jac, void *jac_data); 00260 00261 /* 00262 * ----------------------------------------------------------------- 00263 * Optional outputs from the CVSDIRECT linear solver 00264 * ----------------------------------------------------------------- 00265 * 00266 * CVDlsGetWorkSpace returns the real and integer workspace used 00267 * by the direct linear solver. 00268 * CVDlsGetNumJacEvals returns the number of calls made to the 00269 * Jacobian evaluation routine jac. 00270 * CVDlsGetNumRhsEvals returns the number of calls to the user 00271 * f routine due to finite difference Jacobian 00272 * evaluation. 00273 * CVDlsGetLastFlag returns the last error flag set by any of 00274 * the CVSDIRECT interface functions. 00275 * 00276 * The return value of CVDlsGet* is one of: 00277 * CVDIRECT_SUCCESS if successful 00278 * CVDIRECT_MEM_NULL if the CVODES memory was NULL 00279 * CVDIRECT_LMEM_NULL if the linear solver memory was NULL 00280 * ----------------------------------------------------------------- 00281 */ 00282 00283 SUNDIALS_EXPORT int CVDlsGetWorkSpace(void *cvode_mem, long int *lenrwLS, long int *leniwLS); 00284 SUNDIALS_EXPORT int CVDlsGetNumJacEvals(void *cvode_mem, long int *njevals); 00285 SUNDIALS_EXPORT int CVDlsGetNumRhsEvals(void *cvode_mem, long int *nfevalsLS); 00286 SUNDIALS_EXPORT int CVDlsGetLastFlag(void *cvode_mem, int *flag); 00287 00288 /* 00289 * ----------------------------------------------------------------- 00290 * The following function returns the name of the constant 00291 * associated with a CVSDIRECT return flag 00292 * ----------------------------------------------------------------- 00293 */ 00294 00295 SUNDIALS_EXPORT char *CVDlsGetReturnFlagName(int flag); 00296 00297 /* 00298 * ================================================================= 00299 * PART II: B A C K W A R D P R O B L E M S 00300 * ================================================================= 00301 */ 00302 00303 /* 00304 * ----------------------------------------------------------------- 00305 * FUNCTION TYPES 00306 * ----------------------------------------------------------------- 00307 */ 00308 00309 /* 00310 * ----------------------------------------------------------------- 00311 * Type: CVDlsDenseJacFnB 00312 * ----------------------------------------------------------------- 00313 * A dense Jacobian approximation function jacB for the adjoint 00314 * (backward) problem must have the prototype given below. 00315 * ----------------------------------------------------------------- 00316 */ 00317 00318 typedef int (*CVDlsDenseJacFnB)(int nB, realtype t, 00319 N_Vector y, 00320 N_Vector yB, N_Vector fyB, 00321 DlsMat JB, void *jac_dataB, 00322 N_Vector tmp1B, N_Vector tmp2B, N_Vector tmp3B); 00323 00324 /* 00325 * ----------------------------------------------------------------- 00326 * Type : CVDlsBandJacFnB 00327 * ----------------------------------------------------------------- 00328 * A band Jacobian approximation function jacB for the adjoint 00329 * (backward) problem must have the prototype given below. 00330 * ----------------------------------------------------------------- 00331 */ 00332 00333 typedef int (*CVDlsBandJacFnB)(int nB, int mupperB, int mlowerB, 00334 realtype t, 00335 N_Vector y, 00336 N_Vector yB, N_Vector fyB, 00337 DlsMat JB, void *jac_dataB, 00338 N_Vector tmp1B, N_Vector tmp2B, N_Vector tmp3B); 00339 00340 /* 00341 * ----------------------------------------------------------------- 00342 * EXPORTED FUNCTIONS 00343 * ----------------------------------------------------------------- 00344 */ 00345 00346 /* 00347 * ----------------------------------------------------------------- 00348 * Functions: CVDlsSetJacFnB 00349 * ----------------------------------------------------------------- 00350 * CVDlsSetJacFn specifies the Jacobian function to be used by a 00351 * CVSDIRECT linear solver for the bacward integration phase. 00352 * ----------------------------------------------------------------- 00353 */ 00354 00355 SUNDIALS_EXPORT int CVDlsSetJacFnB(void *cvadj_mem, void *jacB, void *jac_dataB); 00356 00357 #ifdef __cplusplus 00358 } 00359 #endif 00360 00361 #endif