ida_bbdpre.h

Go to the documentation of this file.
00001 /*
00002  * -----------------------------------------------------------------
00003  * $Revision: 1.3 $
00004  * $Date: 2006/11/29 00:05:06 $
00005  * ----------------------------------------------------------------- 
00006  * Programmer(s): Alan C. Hindmarsh, Radu Serban and
00007  *                Aaron Collier @ LLNL
00008  * -----------------------------------------------------------------
00009  * Copyright (c) 2002, The Regents of the University of California.
00010  * Produced at the Lawrence Livermore National Laboratory.
00011  * All rights reserved.
00012  * For details, see the LICENSE file.
00013  * -----------------------------------------------------------------
00014  * This is the header file for the IDABBDPRE module, for a
00015  * band-block-diagonal preconditioner, i.e. a block-diagonal
00016  * matrix with banded blocks, for use with IDA and
00017  * IDASPGMR/IDASPBCG/IDASPTFQMR.
00018  *
00019  * Summary:
00020  *
00021  * These routines provide a preconditioner matrix that is
00022  * block-diagonal with banded blocks. The blocking corresponds
00023  * to the distribution of the dependent variable vector y among
00024  * the processors. Each preconditioner block is generated from
00025  * the Jacobian of the local part (on the current processor) of a
00026  * given function G(t,y,y') approximating F(t,y,y'). The blocks
00027  * are generated by a difference quotient scheme on each processor
00028  * independently. This scheme utilizes an assumed banded structure
00029  * with given half-bandwidths, mudq and mldq. However, the banded
00030  * Jacobian block kept by the scheme has half-bandwiths mukeep and
00031  * mlkeep, which may be smaller.
00032  *
00033  * The user's calling program should have the following form:
00034  *
00035  *   #include <ida/ida_bbdpre.h>
00036  *   #include <nvector_parallel.h>
00037  *   ...
00038  *   void *p_data;
00039  *   ...
00040  *   y0  = N_VNew_Parallel(...);
00041  *   yp0 = N_VNew_Parallel(...);
00042  *   ...
00043  *   ida_mem = IDACreate(...);
00044  *   ier = IDAMalloc(...);
00045  *   ...
00046  *   p_data = IDABBDPrecAlloc(ida_mem, Nlocal, mudq, mldq,
00047  *                            mukeep, mlkeep, dq_rel_yy, Gres, Gcomm);
00048  *   flag = IDABBDSptfqmr(ida_mem, maxl, p_data);
00049  *       -or-
00050  *   flag = IDABBDSpgmr(ida_mem, maxl, p_data);
00051  *       -or-
00052  *   flag = IDABBDSpbcg(ida_mem, maxl, p_data);
00053  *   ...
00054  *   ier = IDASolve(...);
00055  *   ...
00056  *   IDABBDFree(&p_data);
00057  *   ...
00058  *   IDAFree(...);
00059  *
00060  *   N_VDestroy(y0);
00061  *   N_VDestroy(yp0);
00062  *
00063  * The user-supplied routines required are:
00064  *
00065  *   res  is the function F(t,y,y') defining the DAE system to
00066  *   be solved:  F(t,y,y') = 0.
00067  *
00068  *   Gres  is the function defining a local approximation
00069  *   G(t,y,y') to F, for the purposes of the preconditioner.
00070  *
00071  *   Gcomm  is the function performing communication needed
00072  *   for Glocal.
00073  *
00074  * Notes:
00075  *
00076  * 1) This header file is included by the user for the definition
00077  *    of the IBBDPrecData type and for needed function prototypes.
00078  *
00079  * 2) The IDABBDPrecAlloc call includes half-bandwidths mudq and
00080  *    mldq to be used in the approximate Jacobian. They need
00081  *    not be the true half-bandwidths of the Jacobian of the
00082  *    local block of G, when smaller values may provide a greater
00083  *    efficiency. Similarly, mukeep and mlkeep, specifying the
00084  *    bandwidth kept for the approximate Jacobian, need not be
00085  *    the true half-bandwidths. Also, mukeep, mlkeep, mudq, and
00086  *    mldq need not be the same on every processor.
00087  *
00088  * 3) The actual name of the user's res function is passed to
00089  *    IDAMalloc, and the names of the user's Gres and Gcomm
00090  *    functions are passed to IDABBDPrecAlloc.        
00091  *
00092  * 4) The pointer to the user-defined data block res_data, which
00093  *    is set through IDASetRdata is also available to the user
00094  *    in glocal and gcomm.
00095  *
00096  * 5) Optional outputs specific to this module are available by
00097  *    way of routines listed below. These include work space sizes
00098  *    and the cumulative number of glocal calls. The costs
00099  *    associated with this module also include nsetups banded LU
00100  *    factorizations, nsetups gcomm calls, and nps banded
00101  *    backsolve calls, where nsetups and nps are integrator
00102  *    optional outputs.
00103  * -----------------------------------------------------------------
00104  */
00105 
00106 #ifndef _IDABBDPRE_H
00107 #define _IDABBDPRE_H
00108 
00109 #ifdef __cplusplus     /* wrapper to enable C++ usage */
00110 extern "C" {
00111 #endif
00112 
00113 #include <sundials/sundials_nvector.h>
00114 
00115 /* IDABBDPRE return values */
00116 
00117 #define IDABBDPRE_SUCCESS       0
00118 #define IDABBDPRE_PDATA_NULL   -11
00119 #define IDABBDPRE_FUNC_UNRECVR -12
00120 
00121 /*
00122  * -----------------------------------------------------------------
00123  * Type : IDABBDLocalFn
00124  * -----------------------------------------------------------------
00125  * The user must supply a function G(t,y,y') which approximates
00126  * the function F for the system F(t,y,y') = 0, and which is
00127  * computed locally (without interprocess communication).
00128  * (The case where G is mathematically identical to F is allowed.)
00129  * The implementation of this function must have type IDABBDLocalFn.
00130  *
00131  * This function takes as input the independent variable value tt,
00132  * the current solution vector yy, the current solution
00133  * derivative vector yp, and a pointer to the user-defined data
00134  * block res_data. It is to compute the local part of G(t,y,y')
00135  * and store it in the vector gval. (Providing memory for yy and
00136  * gval is handled within this preconditioner module.) It is
00137  * expected that this routine will save communicated data in work
00138  * space defined by the user, and made available to the
00139  * preconditioner function for the problem. The res_data
00140  * parameter is the same as that passed by the user to the
00141  * IDASetRdata routine.
00142  *
00143  * An IDABBDLocalFn Gres is to return an int, defined in the same
00144  * way as for the residual function: 0 (success), +1 or -1 (fail).
00145  * -----------------------------------------------------------------
00146  */
00147 
00148 typedef int (*IDABBDLocalFn)(int Nlocal, realtype tt,
00149                  N_Vector yy, N_Vector yp, N_Vector gval,
00150                  void *res_data);
00151 
00152 /*
00153  * -----------------------------------------------------------------
00154  * Type : IDABBDCommFn
00155  * -----------------------------------------------------------------
00156  * The user may supply a function of type IDABBDCommFn which
00157  * performs all interprocess communication necessary to
00158  * evaluate the approximate system function described above.
00159  *
00160  * This function takes as input the solution vectors yy and yp,
00161  * and a pointer to the user-defined data block res_data. The
00162  * res_data parameter is the same as that passed by the user to
00163  * the IDAMalloc routine.
00164  *
00165  * The IDABBDCommFn Gcomm is expected to save communicated data in
00166  * space defined with the structure *res_data.
00167  *
00168  * A IDABBDCommFn Gcomm returns an int value equal to 0 (success),
00169  * > 0 (recoverable error), or < 0 (unrecoverable error).
00170  *
00171  * Each call to the IDABBDCommFn is preceded by a call to the system
00172  * function res with the same vectors yy and yp. Thus the
00173  * IDABBDCommFn gcomm can omit any communications done by res if
00174  * relevant to the evaluation of the local function glocal.
00175  * A NULL communication function can be passed to IDABBDPrecAlloc
00176  * if all necessary communication was done by res.
00177  * -----------------------------------------------------------------
00178  */
00179 
00180 typedef int (*IDABBDCommFn)(int Nlocal, realtype tt,
00181                 N_Vector yy, N_Vector yp,
00182                 void *res_data);
00183 
00184 /*
00185  * -----------------------------------------------------------------
00186  * Function : IDABBDPrecAlloc
00187  * -----------------------------------------------------------------
00188  * IDABBDPrecAlloc allocates and initializes an IBBDPrecData
00189  * structure to be passed to IDASp* (and used by IDABBDPrecSetup
00190  * and IDABBDPrecSol).
00191  *
00192  * The parameters of IDABBDPrecAlloc are as follows:
00193  *
00194  * ida_mem  is a pointer to the memory blockreturned by IDACreate.
00195  *
00196  * Nlocal  is the length of the local block of the vectors yy etc.
00197  *         on the current processor.
00198  *
00199  * mudq, mldq  are the upper and lower half-bandwidths to be used
00200  *             in the computation of the local Jacobian blocks.
00201  *
00202  * mukeep, mlkeep are the upper and lower half-bandwidths to be
00203  *                used in saving the Jacobian elements in the local
00204  *                block of the preconditioner matrix PP.
00205  *
00206  * dq_rel_yy is an optional input. It is the relative increment
00207  *           to be used in the difference quotient routine for
00208  *           Jacobian calculation in the preconditioner. The
00209  *           default is sqrt(unit roundoff), and specified by
00210  *           passing dq_rel_yy = 0.
00211  *
00212  * Gres    is the name of the user-supplied function G(t,y,y')
00213  *         that approximates F and whose local Jacobian blocks
00214  *         are to form the preconditioner.
00215  *
00216  * Gcomm   is the name of the user-defined function that performs
00217  *         necessary interprocess communication for the
00218  *         execution of glocal.
00219  *
00220  * IDABBDPrecAlloc returns the storage allocated (type *void),
00221  * or NULL if the request for storage cannot be satisfied.
00222  * -----------------------------------------------------------------
00223  */
00224 
00225 SUNDIALS_EXPORT void *IDABBDPrecAlloc(void *ida_mem, int Nlocal,
00226                       int mudq, int mldq,
00227                       int mukeep, int mlkeep,
00228                       realtype dq_rel_yy,
00229                       IDABBDLocalFn Gres, IDABBDCommFn Gcomm);
00230 
00231 /*
00232  * -----------------------------------------------------------------
00233  * Function : IDABBDSptfqmr
00234  * -----------------------------------------------------------------
00235  * IDABBDSptfqmr links the IDABBDPRE preconditioner to the IDASPTFQMR
00236  * linear solver. It performs the following actions:
00237  *  1) Calls the IDASPTFQMR specification routine and attaches the
00238  *     IDASPTFQMR linear solver to the IDA solver;
00239  *  2) Sets the preconditioner data structure for IDASPTFQMR
00240  *  3) Sets the preconditioner setup routine for IDASPTFQMR
00241  *  4) Sets the preconditioner solve routine for IDASPTFQMR
00242  *
00243  * Its first 2 arguments are the same as for IDASptfqmr (see
00244  * idasptfqmr.h). The last argument is the pointer to the IDABBDPRE
00245  * memory block returned by IDABBDPrecAlloc. Note that the user need
00246  * not call IDASptfqmr anymore.
00247  *
00248  * Possible return values are:
00249  *    IDASPILS_SUCCESS     if successful
00250  *    IDASPILS_MEM_NULL    if the IDA memory was NULL
00251  *    IDASPILS_MEM_FAIL    if there was a memory allocation failure
00252  *    IDASPILS_ILL_INPUT   if there was illegal input
00253  *    IDABBDPRE_PDATA_NULL if bbd_data was NULL
00254  * -----------------------------------------------------------------
00255  */
00256 
00257 SUNDIALS_EXPORT int IDABBDSptfqmr(void *ida_mem, int maxl, void *bbd_data);
00258 
00259 /*
00260  * -----------------------------------------------------------------
00261  * Function : IDABBDSpbcg
00262  * -----------------------------------------------------------------
00263  * IDABBDSpbcg links the IDABBDPRE preconditioner to the IDASPBCG
00264  * linear solver. It performs the following actions:
00265  *  1) Calls the IDASPBCG specification routine and attaches the
00266  *     IDASPBCG linear solver to the IDA solver;
00267  *  2) Sets the preconditioner data structure for IDASPBCG
00268  *  3) Sets the preconditioner setup routine for IDASPBCG
00269  *  4) Sets the preconditioner solve routine for IDASPBCG
00270  *
00271  * Its first 2 arguments are the same as for IDASpbcg (see
00272  * idaspbcg.h). The last argument is the pointer to the IDABBDPRE
00273  * memory block returned by IDABBDPrecAlloc. Note that the user need
00274  * not call IDASpbcg anymore.
00275  *
00276  * Possible return values are:
00277  *    IDASPILS_SUCCESS      if successful
00278  *    IDASPILS_MEM_NULL     if the IDA memory was NULL
00279  *    IDASPILS_MEM_FAIL     if there was a memory allocation failure
00280  *    IDASPILS_ILL_INPUT    if there was illegal input
00281  *    IDABBDPRE_PDATA_NULL  if bbd_data was NULL
00282  * -----------------------------------------------------------------
00283  */
00284 
00285 SUNDIALS_EXPORT int IDABBDSpbcg(void *ida_mem, int maxl, void *bbd_data);
00286 
00287 /*
00288  * -----------------------------------------------------------------
00289  * Function : IDABBDSpgmr
00290  * -----------------------------------------------------------------
00291  * IDABBDSpgmr links the IDABBDPRE preconditioner to the IDASPGMR
00292  * linear solver. It performs the following actions:
00293  *  1) Calls the IDASPGMR specification routine and attaches the
00294  *     IDASPGMR linear solver to the IDA solver;
00295  *  2) Sets the preconditioner data structure for IDASPGMR
00296  *  3) Sets the preconditioner setup routine for IDASPGMR
00297  *  4) Sets the preconditioner solve routine for IDASPGMR
00298  *
00299  * Its first 2 arguments are the same as for IDASpgmr (see
00300  * idaspgmr.h). The last argument is the pointer to the IDABBDPRE
00301  * memory block returned by IDABBDPrecAlloc. Note that the user need
00302  * not call IDASpgmr anymore.
00303  *
00304  * Possible return values are:
00305  *    IDASPILS_SUCCESS      if successful
00306  *    IDASPILS_MEM_NULL     if the ida memory was NULL
00307  *    IDASPILS_MEM_FAIL     if there was a memory allocation failure
00308  *    IDASPILS_ILL_INPUT    if there was illegal input
00309  *    IDABBDPRE_PDATA_NULL  if bbd_data was NULL
00310  * -----------------------------------------------------------------
00311  */
00312 
00313 SUNDIALS_EXPORT int IDABBDSpgmr(void *ida_mem, int maxl, void *bbd_data);
00314 
00315 /*
00316  * -----------------------------------------------------------------
00317  * Function : IDABBDPrecReInit
00318  * -----------------------------------------------------------------
00319  * IDABBDPrecReInit reinitializes the IDABBDPRE module when
00320  * solving a sequence of problems of the same size with
00321  * IDASPGMR/IDABBDPRE, IDASPBCG/IDABBDPRE, or IDASPTFQMR/IDABBDPRE
00322  * provided there is no change in Nlocal, mukeep, or mlkeep. After
00323  * solving one problem, and after calling IDAReInit to reinitialize
00324  * the integrator for a subsequent problem, call IDABBDPrecReInit.
00325  *
00326  * The first argument to IDABBDPrecReInit must be the pointer
00327  * bbd_data that was returned by IDABBDPrecAlloc. All other
00328  * arguments have the same names and meanings as those of
00329  * IDABBDPrecAlloc.
00330  *
00331  * The return value of IDABBDPrecReInit is IDABBDPRE_SUCCESS, indicating
00332  * success, or IDABBDPRE_PDATA_NULL if bbd_data was NULL.
00333  * -----------------------------------------------------------------
00334  */
00335 
00336 SUNDIALS_EXPORT int IDABBDPrecReInit(void *bbd_data,
00337                      int mudq, int mldq,
00338                      realtype dq_rel_yy,
00339                      IDABBDLocalFn Gres, IDABBDCommFn Gcomm);
00340 
00341 /*
00342  * -----------------------------------------------------------------
00343  * Function : IDABBDPrecFree
00344  * -----------------------------------------------------------------
00345  * IDABBDPrecFree frees the memory block bbd_data allocated by the
00346  * call to IDABBDPrecAlloc.
00347  * -----------------------------------------------------------------
00348  */
00349 
00350 SUNDIALS_EXPORT void IDABBDPrecFree(void **bbd_data);
00351 
00352 /*
00353  * -----------------------------------------------------------------
00354  * Optional outputs for IDABBDPRE
00355  * -----------------------------------------------------------------
00356  * IDABBDPrecGetWorkSpace returns the real and integer work space
00357  *                        for IDABBDPRE.
00358  * IDABBDPrecGetNumGfnEvals returns the number of calls to the
00359  *                          user Gres function.
00360  * 
00361  * The return value of IDABBDPrecGet* is one of:
00362  *    IDABBDPRE_SUCCESS    if successful
00363  *    IDABBDPRE_PDATA_NULL if the bbd_data memory was NULL
00364  * -----------------------------------------------------------------
00365  */
00366 
00367 SUNDIALS_EXPORT int IDABBDPrecGetWorkSpace(void *bbd_data, 
00368                                            long int *lenrwBBDP, long int *leniwBBDP);
00369 SUNDIALS_EXPORT int IDABBDPrecGetNumGfnEvals(void *bbd_data, long int *ngevalsBBDP);
00370 
00371 /*
00372  * -----------------------------------------------------------------
00373  * The following function returns the name of the constant 
00374  * associated with an IDABBDPRE return flag
00375  * -----------------------------------------------------------------
00376  */
00377   
00378 SUNDIALS_EXPORT char *IDABBDPrecGetReturnFlagName(int flag);
00379 
00380 
00381 #ifdef __cplusplus
00382 }
00383 #endif
00384 
00385 #endif

Generated on Fri Sep 26 07:44:10 2008 for SimTKcore by  doxygen 1.5.6