00001 /* 00002 * ----------------------------------------------------------------- 00003 * $Revision: 1.2 $ 00004 * $Date: 2006/11/29 00:05:07 $ 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 KINSOL. 00014 * ----------------------------------------------------------------- 00015 */ 00016 00017 #ifndef _KINDIRECT_H 00018 #define _KINDIRECT_H 00019 00020 #ifdef __cplusplus /* wrapper to enable C++ usage */ 00021 extern "C" { 00022 #endif 00023 00024 #include <sundials/sundials_direct.h> 00025 #include <sundials/sundials_nvector.h> 00026 00027 /* 00028 * ================================================================= 00029 * K I N D I R E C T C O N S T A N T S 00030 * ================================================================= 00031 */ 00032 00033 /* 00034 * ----------------------------------------------------------------- 00035 * KINDIRECT return values 00036 * ----------------------------------------------------------------- 00037 */ 00038 00039 #define KINDIRECT_SUCCESS 0 00040 #define KINDIRECT_MEM_NULL -1 00041 #define KINDIRECT_LMEM_NULL -2 00042 #define KINDIRECT_ILL_INPUT -3 00043 #define KINDIRECT_MEM_FAIL -4 00044 00045 /* Additional last_flag values */ 00046 00047 #define KINDIRECT_JACFUNC_UNRECVR -5 00048 #define KINDIRECT_JACFUNC_RECVR -6 00049 00050 /* 00051 * ================================================================= 00052 * F U N C T I O N T Y P E S 00053 * ================================================================= 00054 */ 00055 00056 /* 00057 * ----------------------------------------------------------------- 00058 * Type: KINDlsDenseJacFn 00059 * ----------------------------------------------------------------- 00060 * 00061 * A dense Jacobian approximation function Jac must be of type 00062 * KINDlsDenseJacFn. Its parameters are: 00063 * 00064 * N - problem size. 00065 * 00066 * u - current iterate (unscaled) [input] 00067 * 00068 * fu - vector (type N_Vector) containing result of nonlinear 00069 * system function evaluated at current iterate: 00070 * fu = F(u) [input] 00071 * 00072 * J - dense matrix (of type DlsMat) that will be loaded 00073 * by a KINDlsDenseJacFn with an approximation to the 00074 * Jacobian matrix J = (dF_i/dy_j). 00075 * 00076 * jac_data - pointer to user data - the same as the jac_data 00077 * parameter passed to KINDlsSetJacFn. 00078 * 00079 * tmp1, tmp2 - available scratch vectors (volatile storage) 00080 * 00081 * A KINDlsDenseJacFn should return 0 if successful, a positive 00082 * value if a recoverable error occurred, and a negative value if 00083 * an unrecoverable error occurred. 00084 * 00085 * ----------------------------------------------------------------- 00086 * 00087 * NOTE: The following are two efficient ways to load a dense Jac: 00088 * (1) (with macros - no explicit data structure references) 00089 * for (j=0; j < Neq; j++) { 00090 * col_j = DENSE_COL(Jac,j); 00091 * for (i=0; i < Neq; i++) { 00092 * generate J_ij = the (i,j)th Jacobian element 00093 * col_j[i] = J_ij; 00094 * } 00095 * } 00096 * (2) (without macros - explicit data structure references) 00097 * for (j=0; j < Neq; j++) { 00098 * col_j = (Jac->data)[j]; 00099 * for (i=0; i < Neq; i++) { 00100 * generate J_ij = the (i,j)th Jacobian element 00101 * col_j[i] = J_ij; 00102 * } 00103 * } 00104 * A third way, using the DENSE_ELEM(A,i,j) macro, is much less 00105 * efficient in general. It is only appropriate for use in small 00106 * problems in which efficiency of access is NOT a major concern. 00107 * 00108 * ----------------------------------------------------------------- 00109 */ 00110 00111 00112 typedef int (*KINDlsDenseJacFn)(int N, 00113 N_Vector u, N_Vector fu, 00114 DlsMat J, void *jac_data, 00115 N_Vector tmp1, N_Vector tmp2); 00116 00117 /* 00118 * ----------------------------------------------------------------- 00119 * Type: KINDlsBandJacFn 00120 * ----------------------------------------------------------------- 00121 * 00122 * A band Jacobian approximation function Jac must have the 00123 * prototype given below. Its parameters are: 00124 * 00125 * N is the problem size 00126 * 00127 * mupper is the upper half-bandwidth of the approximate banded 00128 * Jacobian. This parameter is the same as the mupper parameter 00129 * passed by the user to the linear solver initialization function. 00130 * 00131 * mlower is the lower half-bandwidth of the approximate banded 00132 * Jacobian. This parameter is the same as the mlower parameter 00133 * passed by the user to the linear solver initialization function. 00134 * 00135 * u - current iterate (unscaled) [input] 00136 * 00137 * fu - vector (type N_Vector) containing result of nonlinear 00138 * system function evaluated at current iterate: 00139 * fu = F(uu) [input] 00140 * 00141 * J - band matrix (of type DlsMat) that will be loaded by a 00142 * KINDlsBandJacFn with an approximation to the Jacobian 00143 * matrix Jac = (dF_i/dy_j). 00144 * 00145 * jac_data - pointer to user data - the same as the jac_data 00146 * parameter passed to KINBandSetJacFn. 00147 * 00148 * tmp1, tmp2 - available scratch vectors (volatile storage) 00149 * 00150 * A KINDlsBandJacFn should return 0 if successful, a positive value 00151 * if a recoverable error occurred, and a negative value if an 00152 * unrecoverable error occurred. 00153 * 00154 * ----------------------------------------------------------------- 00155 * 00156 * NOTE. Three efficient ways to load J are: 00157 * 00158 * (1) (with macros - no explicit data structure references) 00159 * for (j=0; j < n; j++) { 00160 * col_j = BAND_COL(Jac,j); 00161 * for (i=j-mupper; i <= j+mlower; i++) { 00162 * generate J_ij = the (i,j)th Jacobian element 00163 * BAND_COL_ELEM(col_j,i,j) = J_ij; 00164 * } 00165 * } 00166 * 00167 * (2) (with BAND_COL macro, but without BAND_COL_ELEM macro) 00168 * for (j=0; j < n; j++) { 00169 * col_j = BAND_COL(Jac,j); 00170 * for (k=-mupper; k <= mlower; k++) { 00171 * generate J_ij = the (i,j)th Jacobian element, i=j+k 00172 * col_j[k] = J_ij; 00173 * } 00174 * } 00175 * 00176 * (3) (without macros - explicit data structure references) 00177 * offset = Jac->smu; 00178 * for (j=0; j < n; j++) { 00179 * col_j = ((Jac->data)[j])+offset; 00180 * for (k=-mupper; k <= mlower; k++) { 00181 * generate J_ij = the (i,j)th Jacobian element, i=j+k 00182 * col_j[k] = J_ij; 00183 * } 00184 * } 00185 * Caution: Jac->smu is generally NOT the same as mupper. 00186 * 00187 * The BAND_ELEM(A,i,j) macro is appropriate for use in small 00188 * problems in which efficiency of access is NOT a major concern. 00189 * 00190 * ----------------------------------------------------------------- 00191 */ 00192 00193 typedef int (*KINDlsBandJacFn)(int N, int mupper, int mlower, 00194 N_Vector u, N_Vector fu, 00195 DlsMat J, void *jac_data, 00196 N_Vector tmp1, N_Vector tmp2); 00197 00198 /* 00199 * ================================================================= 00200 * E X P O R T E D F U N C T I O N S 00201 * ================================================================= 00202 */ 00203 00204 /* 00205 * ----------------------------------------------------------------- 00206 * Optional inputs to the KINDIRECT linear solver 00207 * ----------------------------------------------------------------- 00208 * 00209 * KINDlsSetJacFn specifies the Jacobian approximation routine to be 00210 * used. When using dense Jacobians, a user-supplied jac routine must 00211 * be of type KINDlsDenseJacFn. When using banded Jacobians, a 00212 * user-supplied jac routine must be of type KINDlsBandJacFn. 00213 * By default, a difference quotient approximation, supplied with this 00214 * solver is used. 00215 * KINDlsSetJacFn also specifies a pointer to user data which is 00216 * passed to the user's jac routine every time it is called. 00217 * 00218 * The return value of KINDlsSetJacFn is one of: 00219 * KINDIRECT_SUCCESS if successful 00220 * KINDIRECT_MEM_NULL if the KINSOL memory was NULL 00221 * KINDIRECT_LMEM_NULL if the linear solver memory was NULL 00222 * ----------------------------------------------------------------- 00223 */ 00224 00225 SUNDIALS_EXPORT int KINDlsSetJacFn(void *kinmem, void *jac, void *jac_data); 00226 00227 /* 00228 * ----------------------------------------------------------------- 00229 * Optional outputs from a KINDIRECT linear solver 00230 * ----------------------------------------------------------------- 00231 * 00232 * KINDlsGetWorkSpace returns the real and integer workspace used 00233 * by the KINDIRECT linear solver. 00234 * KINDlsGetNumJacEvals returns the number of calls made to the 00235 * Jacobian evaluation routine. 00236 * KINDlsGetNumFuncEvals returns the number of calls to the user's F 00237 * routine due to finite difference Jacobian 00238 * evaluation. 00239 * KINDlsGetLastFlag returns the last error flag set by any of 00240 * the KINDIRECT interface functions. 00241 * KINDlsGetReturnFlagName returns the name of the constant 00242 * associated with a KINDIRECT return flag 00243 * 00244 * The return value of KINDlsGet* is one of: 00245 * KINDIRECT_SUCCESS if successful 00246 * KINDIRECT_MEM_NULL if the KINSOL memory was NULL 00247 * KINDIRECT_LMEM_NULL if the linear solver memory was NULL 00248 * ----------------------------------------------------------------- 00249 */ 00250 00251 SUNDIALS_EXPORT int KINDlsGetWorkSpace(void *kinmem, long int *lenrwB, long int *leniwB); 00252 SUNDIALS_EXPORT int KINDlsGetNumJacEvals(void *kinmem, long int *njevalsB); 00253 SUNDIALS_EXPORT int KINDlsGetNumFuncEvals(void *kinmem, long int *nfevalsB); 00254 SUNDIALS_EXPORT int KINDlsGetLastFlag(void *kinmem, int *flag); 00255 SUNDIALS_EXPORT char *KINDlsGetReturnFlagName(int flag); 00256 00257 #ifdef __cplusplus 00258 } 00259 #endif 00260 00261 #endif