00001 /* 00002 * ----------------------------------------------------------------- 00003 * $Revision: 1.3 $ 00004 * $Date: 2006/11/29 00:05:07 $ 00005 * ----------------------------------------------------------------- 00006 * Programmer(s): Alan Hindmarsh, Radu Serban, and Aaron Collier @ LLNL 00007 * ----------------------------------------------------------------- 00008 * Copyright (c) 2002, 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 * This is the header file for the KINBBDPRE module, for a 00014 * band-block-diagonal preconditioner, i.e. a block-diagonal 00015 * matrix with banded blocks, for use with KINSol, KINSp*, 00016 * and the parallel implementaion of the NVECTOR module. 00017 * 00018 * Summary: 00019 * 00020 * These routines provide a preconditioner matrix for KINSol that 00021 * is block-diagonal with banded blocks. The blocking corresponds 00022 * to the distribution of the dependent variable vector u amongst 00023 * the processes. Each preconditioner block is generated from 00024 * the Jacobian of the local part (associated with the current 00025 * process) of a given function g(u) approximating f(u). The blocks 00026 * are generated by each process via a difference quotient scheme, 00027 * utilizing the assumed banded structure with given half-bandwidths, 00028 * mudq and mldq. However, the banded Jacobian block kept by the 00029 * scheme has half-bandwidths mukeep and mlkeep, which may be smaller. 00030 * 00031 * The user's calling program should have the following form: 00032 * 00033 * #include <sundials/sundials_types.h> 00034 * #include <sundials/sundials_math.h> 00035 * #include <sundials/sundials_iterative.h> 00036 * #include <nvector_parallel.h> 00037 * #include <kinsol.h> 00038 * #include <kinsol/kinsol_bbdpre.h> 00039 * ... 00040 * void *p_data; 00041 * ... 00042 * MPI_Init(&argc,&argv); 00043 * ... 00044 * tmpl = N_VNew_Parallel(...); 00045 * ... 00046 * kin_mem = KINCreate(); 00047 * KINMalloc(kin_mem,...,tmpl); 00048 * ... 00049 * p_data = KINBBDPrecAlloc(kin_mem,...); 00050 * ... 00051 * KINBBDSptfqmr(kin_mem,...,p_data); 00052 * -or- 00053 * KINBBDSpbcg(kin_mem,...,p_data); 00054 * -or- 00055 * KINBBDSpgmr(kin_mem,...,p_data); 00056 * ... 00057 * KINSol(kin_mem,...); 00058 * ... 00059 * KINBBDPrecFree(&p_data); 00060 * ... 00061 * KINFree(kin_mem); 00062 * ... 00063 * N_VDestroy_Parallel(tmpl); 00064 * ... 00065 * MPI_Finalize(); 00066 * 00067 * The user-supplied routines required are: 00068 * 00069 * func the function f(u) defining the system to be solved: 00070 * f(u) = 0 00071 * 00072 * glocal the function defining the approximation g(u) to f(u) 00073 * 00074 * gcomm the function to do necessary communication for glocal 00075 * 00076 * Notes: 00077 * 00078 * 1) This header file (kinsol_bbdpre.h) is included by the user for 00079 * the definition of the KBBDData data type and for needed 00080 * function prototypes. 00081 * 00082 * 2) The KINBBDPrecAlloc call includes half-bandwiths mudq and mldq 00083 * to be used in the difference quotient calculation of the 00084 * approximate Jacobian. They need not be the true half-bandwidths 00085 * of the Jacobian of the local block of g, when smaller values may 00086 * provide greater efficiency. Also, the half-bandwidths mukeep and 00087 * mlkeep of the retained banded approximate Jacobian block may be 00088 * even smaller, to furhter reduce storage and computational costs. 00089 * For all four half-bandwidths, the values need not be the same 00090 * for every process. 00091 * 00092 * 3) The actual name of the user's f function is passed to 00093 * KINMalloc, and the names of the user's glocal and gcomm 00094 * functions are passed to KINBBDPrecAlloc. 00095 * 00096 * 4) Optional outputs specific to this module are available by 00097 * way of the functions listed below. These include work space 00098 * sizes and the cumulative number of glocal calls. 00099 * ----------------------------------------------------------------- 00100 */ 00101 00102 #ifndef _KINBBDPRE_H 00103 #define _KINBBDPRE_H 00104 00105 #ifdef __cplusplus /* wrapper to enable C++ usage */ 00106 extern "C" { 00107 #endif 00108 00109 #include <sundials/sundials_nvector.h> 00110 00111 /* KINBBDPRE return values */ 00112 00113 #define KINBBDPRE_SUCCESS 0 00114 #define KINBBDPRE_PDATA_NULL -11 00115 #define KINBBDPRE_FUNC_UNRECVR -12 00116 /* 00117 * ----------------------------------------------------------------- 00118 * Type : KINCommFn 00119 * ----------------------------------------------------------------- 00120 * The user must supply a function of type KINCommFn which 00121 * performs all inter-process communication necessary to 00122 * evaluate the approximate system function described above. 00123 * 00124 * This function takes as input the local vector size Nlocal, 00125 * the solution vector u, and a pointer to the user-defined 00126 * data block f_data. 00127 * 00128 * The KINCommFn gcomm is expected to save communicated data in 00129 * space defined with the structure *f_data. 00130 * 00131 * Each call to the KINCommFn is preceded by a call to the system 00132 * function func at the current iterate uu. Thus functions of the 00133 * type KINCommFn can omit any communications done by f (func) if 00134 * relevant to the evaluation of the KINLocalFn function. If all 00135 * necessary communication was done in func, the user can pass 00136 * NULL for gcomm in the call to KINBBDPrecAlloc (see below). 00137 * 00138 * A KINCommFn function should return 0 if successful or 00139 * a non-zero value if an error occured. 00140 * ----------------------------------------------------------------- 00141 */ 00142 00143 typedef int (*KINCommFn)(int Nlocal, N_Vector u, void *f_data); 00144 00145 /* 00146 * ----------------------------------------------------------------- 00147 * Type : KINLocalFn 00148 * ----------------------------------------------------------------- 00149 * The user must supply a function g(u) which approximates the 00150 * function f for the system f(u) = 0, and which is computed 00151 * locally (without inter-process communication). Note: The case 00152 * where g is mathematically identical to f is allowed. 00153 * 00154 * The implementation of this function must have type KINLocalFn 00155 * and take as input the local vector size Nlocal, the local 00156 * solution vector uu, the returned local g values vector, and a 00157 * pointer to the user-defined data block f_data. It is to 00158 * compute the local part of g(u) and store the result in the 00159 * vector gval. (Note: Memory for uu and gval is handled within the 00160 * preconditioner module.) It is expected that this routine will 00161 * save communicated data in work space defined by the user and 00162 * made available to the preconditioner function for the problem. 00163 * 00164 * A KINLocalFn function should return 0 if successful or 00165 * a non-zero value if an error occured. 00166 * ----------------------------------------------------------------- 00167 */ 00168 00169 typedef int (*KINLocalFn)(int Nlocal, N_Vector uu, 00170 N_Vector gval, void *f_data); 00171 00172 /* 00173 * ----------------------------------------------------------------- 00174 * Function : KINBBDPrecAlloc 00175 * ----------------------------------------------------------------- 00176 * KINBBDPrecAlloc allocates and initializes a KBBDData data 00177 * structure to be passed to KINSpgmr/KINSpbcg (and then used by 00178 * KINBBDPrecSetup and KINBBDPrecSolve). 00179 * 00180 * The parameters of KINBBDPrecAlloc are as follows: 00181 * 00182 * kinmem is a pointer to the KINSol memory block. 00183 * 00184 * Nlocal is the length of the local block of the vectors 00185 * on the current process. 00186 * 00187 * mudq, mldq are the upper and lower half-bandwidths to be used 00188 * in the computation of the local Jacobian blocks. 00189 * 00190 * mukeep, mlkeep are the upper and lower half-bandwidths of the 00191 * retained banded approximation to the local 00192 * Jacobian block. 00193 * 00194 * dq_rel_uu is the relative error to be used in the difference 00195 * quotient Jacobian calculation in the preconditioner. 00196 * The default is sqrt(unit roundoff), obtained by 00197 * passing 0. 00198 * 00199 * gloc is the name of the user-supplied function g(u) that 00200 * approximates f and whose local Jacobian blocks are 00201 * to form the preconditioner. 00202 * 00203 * gcomm is the name of the user-defined function that performs 00204 * necessary inter-process communication for the 00205 * execution of gloc. 00206 * 00207 * Note: KINBBDPrecAlloc returns a pointer to the storage allocated, 00208 * or NULL if the request for storage cannot be satisfied. 00209 * ----------------------------------------------------------------- 00210 */ 00211 00212 SUNDIALS_EXPORT void *KINBBDPrecAlloc(void *kinmem, int Nlocal, 00213 int mudq, int mldq, 00214 int mukeep, int mlkeep, 00215 realtype dq_rel_uu, 00216 KINLocalFn gloc, KINCommFn gcomm); 00217 00218 /* 00219 * ----------------------------------------------------------------- 00220 * Function : KINBBDSptfqmr 00221 * ----------------------------------------------------------------- 00222 * KINBBDSptfqmr links the KINBBDPRE preconditioner to the KINSPTFQMR 00223 * linear solver. It performs the following actions: 00224 * 1) Calls the KINSPTFQMR specification routine and attaches the 00225 * KINSPTFQMR linear solver to the KINSOL solver; 00226 * 2) Sets the preconditioner data structure for KINSPTFQMR 00227 * 3) Sets the preconditioner setup routine for KINSPTFQMR 00228 * 4) Sets the preconditioner solve routine for KINSPTFQMR 00229 * 00230 * Its first 2 arguments are the same as for KINSptfqmr (see 00231 * kinsptfqmr.h). The last argument is the pointer to the KBBDPRE 00232 * memory block returned by KINBBDPrecAlloc. 00233 * 00234 * Possible return values are: 00235 * (from kinsptfqmr.h) KINSPTFQMR_SUCCESS 00236 * KINSPTFQMR_MEM_NULL 00237 * KINSPTFQMR_LMEM_NULL 00238 * KINSPTFQMR_MEM_FAIL 00239 * KINSPTFQMR_ILL_INPUT 00240 * 00241 * Additionaly, if KINBBDPrecAlloc was not previously called, 00242 * KINBBDSptfqmr returns KINBBDPRE_PDATA_NULL. 00243 * ----------------------------------------------------------------- 00244 */ 00245 00246 SUNDIALS_EXPORT int KINBBDSptfqmr(void *kinmem, int maxl, void *p_data); 00247 00248 /* 00249 * ----------------------------------------------------------------- 00250 * Function : KINBBDSpbcg 00251 * ----------------------------------------------------------------- 00252 * KINBBDSpbcg links the KINBBDPRE preconditioner to the KINSPBCG 00253 * linear solver. It performs the following actions: 00254 * 1) Calls the KINSPBCG specification routine and attaches the 00255 * KINSPBCG linear solver to the KINSOL solver; 00256 * 2) Sets the preconditioner data structure for KINSPBCG 00257 * 3) Sets the preconditioner setup routine for KINSPBCG 00258 * 4) Sets the preconditioner solve routine for KINSPBCG 00259 * 00260 * Its first 2 arguments are the same as for KINSpbcg (see 00261 * kinspbcg.h). The last argument is the pointer to the KBBDPRE 00262 * memory block returned by KINBBDPrecAlloc. 00263 * 00264 * Possible return values are: 00265 * (from kinspbcg.h) KINSPBCG_SUCCESS 00266 * KINSPBCG_MEM_NULL 00267 * KINSPBCG_LMEM_NULL 00268 * KINSPBCG_MEM_FAIL 00269 * KINSPBCG_ILL_INPUT 00270 * 00271 * Additionaly, if KINBBDPrecAlloc was not previously called, 00272 * KINBBDSpbcg returns KINBBDPRE_PDATA_NULL. 00273 * ----------------------------------------------------------------- 00274 */ 00275 00276 SUNDIALS_EXPORT int KINBBDSpbcg(void *kinmem, int maxl, void *p_data); 00277 00278 /* 00279 * ----------------------------------------------------------------- 00280 * Function : KINBBDSpgmr 00281 * ----------------------------------------------------------------- 00282 * KINBBDSpgmr links the KINBBDPRE preconditioner to the KINSPGMR 00283 * linear solver. It performs the following actions: 00284 * 1) Calls the KINSPGMR specification routine and attaches the 00285 * KINSPGMR linear solver to the KINSOL solver; 00286 * 2) Sets the preconditioner data structure for KINSPGMR 00287 * 3) Sets the preconditioner setup routine for KINSPGMR 00288 * 4) Sets the preconditioner solve routine for KINSPGMR 00289 * 00290 * Its first 2 arguments are the same as for KINSpgmr (see 00291 * kinspgmr.h). The last argument is the pointer to the KBBDPRE 00292 * memory block returned by KINBBDPrecAlloc. 00293 * 00294 * Possible return values are: 00295 * (from kinspgmr.h) KINSPGMR_SUCCESS 00296 * KINSPGMR_MEM_NULL 00297 * KINSPGMR_LMEM_NULL 00298 * KINSPGMR_MEM_FAIL 00299 * KINSPGMR_ILL_INPUT 00300 * 00301 * Additionaly, if KINBBDPrecAlloc was not previously called, 00302 * KINBBDSpgmr returns KINBBDPRE_PDATA_NULL. 00303 * ----------------------------------------------------------------- 00304 */ 00305 00306 SUNDIALS_EXPORT int KINBBDSpgmr(void *kinmem, int maxl, void *p_data); 00307 00308 /* 00309 * ----------------------------------------------------------------- 00310 * Function : KINBBDPrecFree 00311 * ----------------------------------------------------------------- 00312 * KINBBDPrecFree frees the memory block p_data allocated by the 00313 * call to KINBBDPrecAlloc. 00314 * ----------------------------------------------------------------- 00315 */ 00316 00317 SUNDIALS_EXPORT void KINBBDPrecFree(void **p_data); 00318 00319 /* 00320 * ----------------------------------------------------------------- 00321 * Function : KINBBDPrecGet* 00322 * 00323 * The return value of KINBBDPrecGet* is one of: 00324 * KINBBDPRE_SUCCESS if successful 00325 * KINBBDPRE_PDATA_NULL if the p_data memory was NULL 00326 * ----------------------------------------------------------------- 00327 */ 00328 00329 SUNDIALS_EXPORT int KINBBDPrecGetWorkSpace(void *p_data, long int *lenrwBBDP, long int *leniwBBDP); 00330 SUNDIALS_EXPORT int KINBBDPrecGetNumGfnEvals(void *p_data, long int *ngevalsBBDP); 00331 00332 /* 00333 * ----------------------------------------------------------------- 00334 * The following function returns the name of the constant 00335 * associated with a KINBBDPRE return flag 00336 * ----------------------------------------------------------------- 00337 */ 00338 00339 SUNDIALS_EXPORT char *KINBBDPrecGetReturnFlagName(int flag); 00340 00341 /* prototypes for functions KINBBDPrecSetup and KINBBDPrecSolve */ 00342 00343 SUNDIALS_EXPORT int KINBBDPrecSetup(N_Vector uu, N_Vector uscale, 00344 N_Vector fval, N_Vector fscale, 00345 void *p_data, 00346 N_Vector vtemp1, N_Vector vtemp2); 00347 00348 SUNDIALS_EXPORT int KINBBDPrecSolve(N_Vector uu, N_Vector uscale, 00349 N_Vector fval, N_Vector fscale, 00350 N_Vector vv, void *p_data, 00351 N_Vector vtemp); 00352 00353 #ifdef __cplusplus 00354 } 00355 #endif 00356 00357 #endif