idas_bbdpre.h 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334
  1. /*
  2. * -----------------------------------------------------------------
  3. * $Revision: 4489 $
  4. * $Date: 2015-04-29 17:15:44 -0700 (Wed, 29 Apr 2015) $
  5. * -----------------------------------------------------------------
  6. * Programmer(s): Alan C. Hindmarsh, Radu Serban and
  7. * Aaron Collier @ LLNL
  8. * -----------------------------------------------------------------
  9. * LLNS Copyright Start
  10. * Copyright (c) 2014, Lawrence Livermore National Security
  11. * This work was performed under the auspices of the U.S. Department
  12. * of Energy by Lawrence Livermore National Laboratory in part under
  13. * Contract W-7405-Eng-48 and in part under Contract DE-AC52-07NA27344.
  14. * Produced at the Lawrence Livermore National Laboratory.
  15. * All rights reserved.
  16. * For details, see the LICENSE file.
  17. * LLNS Copyright End
  18. * -----------------------------------------------------------------
  19. * This is the header file for the IDABBDPRE module, for a
  20. * band-block-diagonal preconditioner, i.e. a block-diagonal
  21. * matrix with banded blocks, for use with IDAS and
  22. * IDASpgmr/IDASpbcg/IDASptfqmr.
  23. */
  24. #ifndef _IDASBBDPRE_H
  25. #define _IDASBBDPRE_H
  26. #include <sundials/sundials_nvector.h>
  27. #ifdef __cplusplus /* wrapper to enable C++ usage */
  28. extern "C" {
  29. #endif
  30. /*
  31. * =================================================================
  32. * PART I - forward problems
  33. * =================================================================
  34. */
  35. /*
  36. * -----------------------------------------------------------------
  37. *
  38. * SUMMARY
  39. *
  40. * These routines provide a preconditioner matrix that is
  41. * block-diagonal with banded blocks. The blocking corresponds
  42. * to the distribution of the dependent variable vector y among
  43. * the processors. Each preconditioner block is generated from
  44. * the Jacobian of the local part (on the current processor) of a
  45. * given function G(t,y,y') approximating F(t,y,y'). The blocks
  46. * are generated by a difference quotient scheme on each processor
  47. * independently. This scheme utilizes an assumed banded structure
  48. * with given half-bandwidths, mudq and mldq. However, the banded
  49. * Jacobian block kept by the scheme has half-bandwiths mukeep and
  50. * mlkeep, which may be smaller.
  51. *
  52. * The user's calling program should have the following form:
  53. *
  54. * #include <idas/idas_bbdpre.h>
  55. * #include <nvector_parallel.h>
  56. * ...
  57. * y0 = N_VNew_Parallel(...);
  58. * yp0 = N_VNew_Parallel(...);
  59. * ...
  60. * ida_mem = IDACreate(...);
  61. * ier = IDAInit(...);
  62. * ...
  63. * flag = IDASptfqmr(ida_mem, maxl);
  64. * -or-
  65. * flag = IDASpgmr(ida_mem, maxl);
  66. * -or-
  67. * flag = IDASpbcg(ida_mem, maxl);
  68. * ...
  69. * flag = IDABBDPrecInit(ida_mem, Nlocal, mudq, mldq,
  70. * mukeep, mlkeep, dq_rel_yy, Gres, Gcomm);
  71. * ...
  72. * ier = IDASolve(...);
  73. * ...
  74. * IDAFree(&ida_mem);
  75. *
  76. * N_VDestroy(y0);
  77. * N_VDestroy(yp0);
  78. *
  79. * The user-supplied routines required are:
  80. *
  81. * res is the function F(t,y,y') defining the DAE system to
  82. * be solved: F(t,y,y') = 0.
  83. *
  84. * Gres is the function defining a local approximation
  85. * G(t,y,y') to F, for the purposes of the preconditioner.
  86. *
  87. * Gcomm is the function performing communication needed
  88. * for Glocal.
  89. *
  90. * Notes:
  91. *
  92. * 1) This header file is included by the user for the definition
  93. * of the IBBDPrecData type and for needed function prototypes.
  94. *
  95. * 2) The IDABBDPrecInit call includes half-bandwidths mudq and
  96. * mldq to be used in the approximate Jacobian. They need
  97. * not be the true half-bandwidths of the Jacobian of the
  98. * local block of G, when smaller values may provide a greater
  99. * efficiency. Similarly, mukeep and mlkeep, specifying the
  100. * bandwidth kept for the approximate Jacobian, need not be
  101. * the true half-bandwidths. Also, mukeep, mlkeep, mudq, and
  102. * mldq need not be the same on every processor.
  103. *
  104. * 3) The actual name of the user's res function is passed to
  105. * IDAInit, and the names of the user's Gres and Gcomm
  106. * functions are passed to IDABBDPrecInit.
  107. *
  108. * 4) The pointer to the user-defined data block user_data, which
  109. * is set through IDASetUserData is also available to the user
  110. * in glocal and gcomm.
  111. *
  112. * 5) Optional outputs specific to this module are available by
  113. * way of routines listed below. These include work space sizes
  114. * and the cumulative number of glocal calls. The costs
  115. * associated with this module also include nsetups banded LU
  116. * factorizations, nsetups gcomm calls, and nps banded
  117. * backsolve calls, where nsetups and nps are integrator
  118. * optional outputs.
  119. * -----------------------------------------------------------------
  120. */
  121. /*
  122. * -----------------------------------------------------------------
  123. * Type : IDABBDLocalFn
  124. * -----------------------------------------------------------------
  125. * The user must supply a function G(t,y,y') which approximates
  126. * the function F for the system F(t,y,y') = 0, and which is
  127. * computed locally (without interprocess communication).
  128. * (The case where G is mathematically identical to F is allowed.)
  129. * The implementation of this function must have type IDABBDLocalFn.
  130. *
  131. * This function takes as input the independent variable value tt,
  132. * the current solution vector yy, the current solution
  133. * derivative vector yp, and a pointer to the user-defined data
  134. * block user_data. It is to compute the local part of G(t,y,y')
  135. * and store it in the vector gval. (Providing memory for yy and
  136. * gval is handled within this preconditioner module.) It is
  137. * expected that this routine will save communicated data in work
  138. * space defined by the user, and made available to the
  139. * preconditioner function for the problem. The user_data
  140. * parameter is the same as that passed by the user to the
  141. * IDAMalloc routine.
  142. *
  143. * An IDABBDLocalFn Gres is to return an int, defined in the same
  144. * way as for the residual function: 0 (success), +1 or -1 (fail).
  145. * -----------------------------------------------------------------
  146. */
  147. typedef int (*IDABBDLocalFn)(long int Nlocal, realtype tt,
  148. N_Vector yy, N_Vector yp, N_Vector gval,
  149. void *user_data);
  150. /*
  151. * -----------------------------------------------------------------
  152. * Type : IDABBDCommFn
  153. * -----------------------------------------------------------------
  154. * The user may supply a function of type IDABBDCommFn which
  155. * performs all interprocess communication necessary to
  156. * evaluate the approximate system function described above.
  157. *
  158. * This function takes as input the solution vectors yy and yp,
  159. * and a pointer to the user-defined data block user_data. The
  160. * user_data parameter is the same as that passed by the user to
  161. * the IDASetUserData routine.
  162. *
  163. * The IDABBDCommFn Gcomm is expected to save communicated data in
  164. * space defined with the structure *user_data.
  165. *
  166. * A IDABBDCommFn Gcomm returns an int value equal to 0 (success),
  167. * > 0 (recoverable error), or < 0 (unrecoverable error).
  168. *
  169. * Each call to the IDABBDCommFn is preceded by a call to the system
  170. * function res with the same vectors yy and yp. Thus the
  171. * IDABBDCommFn gcomm can omit any communications done by res if
  172. * relevant to the evaluation of the local function glocal.
  173. * A NULL communication function can be passed to IDABBDPrecInit
  174. * if all necessary communication was done by res.
  175. * -----------------------------------------------------------------
  176. */
  177. typedef int (*IDABBDCommFn)(long int Nlocal, realtype tt,
  178. N_Vector yy, N_Vector yp,
  179. void *user_data);
  180. /*
  181. * -----------------------------------------------------------------
  182. * Function : IDABBDPrecInit
  183. * -----------------------------------------------------------------
  184. * IDABBDPrecInit allocates and initializes the BBD preconditioner.
  185. *
  186. * The parameters of IDABBDPrecInit are as follows:
  187. *
  188. * ida_mem is a pointer to the memory blockreturned by IDACreate.
  189. *
  190. * Nlocal is the length of the local block of the vectors yy etc.
  191. * on the current processor.
  192. *
  193. * mudq, mldq are the upper and lower half-bandwidths to be used
  194. * in the computation of the local Jacobian blocks.
  195. *
  196. * mukeep, mlkeep are the upper and lower half-bandwidths to be
  197. * used in saving the Jacobian elements in the local
  198. * block of the preconditioner matrix PP.
  199. *
  200. * dq_rel_yy is an optional input. It is the relative increment
  201. * to be used in the difference quotient routine for
  202. * Jacobian calculation in the preconditioner. The
  203. * default is sqrt(unit roundoff), and specified by
  204. * passing dq_rel_yy = 0.
  205. *
  206. * Gres is the name of the user-supplied function G(t,y,y')
  207. * that approximates F and whose local Jacobian blocks
  208. * are to form the preconditioner.
  209. *
  210. * Gcomm is the name of the user-defined function that performs
  211. * necessary interprocess communication for the
  212. * execution of glocal.
  213. *
  214. * The return value of IDABBDPrecInit is one of:
  215. * IDASPILS_SUCCESS if no errors occurred
  216. * IDASPILS_MEM_NULL if the integrator memory is NULL
  217. * IDASPILS_LMEM_NULL if the linear solver memory is NULL
  218. * IDASPILS_ILL_INPUT if an input has an illegal value
  219. * IDASPILS_MEM_FAIL if a memory allocation request failed
  220. * -----------------------------------------------------------------
  221. */
  222. SUNDIALS_EXPORT int IDABBDPrecInit(void *ida_mem, long int Nlocal,
  223. long int mudq, long int mldq,
  224. long int mukeep, long int mlkeep,
  225. realtype dq_rel_yy,
  226. IDABBDLocalFn Gres, IDABBDCommFn Gcomm);
  227. /*
  228. * -----------------------------------------------------------------
  229. * Function : IDABBDPrecReInit
  230. * -----------------------------------------------------------------
  231. * IDABBDPrecReInit reinitializes the IDABBDPRE module when
  232. * solving a sequence of problems of the same size with
  233. * IDASPGMR/IDABBDPRE, IDASPBCG/IDABBDPRE, or IDASPTFQMR/IDABBDPRE
  234. * provided there is no change in Nlocal, mukeep, or mlkeep. After
  235. * solving one problem, and after calling IDAReInit to reinitialize
  236. * the integrator for a subsequent problem, call IDABBDPrecReInit.
  237. *
  238. * All arguments have the same names and meanings as those
  239. * of IDABBDPrecInit.
  240. *
  241. * The return value of IDABBDPrecReInit is one of:
  242. * IDASPILS_SUCCESS if no errors occurred
  243. * IDASPILS_MEM_NULL if the integrator memory is NULL
  244. * IDASPILS_LMEM_NULL if the linear solver memory is NULL
  245. * IDASPILS_PMEM_NULL if the preconditioner memory is NULL
  246. * -----------------------------------------------------------------
  247. */
  248. SUNDIALS_EXPORT int IDABBDPrecReInit(void *ida_mem,
  249. long int mudq, long int mldq,
  250. realtype dq_rel_yy);
  251. /*
  252. * -----------------------------------------------------------------
  253. * Optional outputs for IDABBDPRE
  254. * -----------------------------------------------------------------
  255. * IDABBDPrecGetWorkSpace returns the real and integer work space
  256. * for IDABBDPRE.
  257. * IDABBDPrecGetNumGfnEvals returns the number of calls to the
  258. * user Gres function.
  259. *
  260. * The return value of IDABBDPrecGet* is one of:
  261. * IDASPILS_SUCCESS if no errors occurred
  262. * IDASPILS_MEM_NULL if the integrator memory is NULL
  263. * IDASPILS_LMEM_NULL if the linear solver memory is NULL
  264. * IDASPILS_PMEM_NULL if the preconditioner memory is NULL
  265. * -----------------------------------------------------------------
  266. */
  267. SUNDIALS_EXPORT int IDABBDPrecGetWorkSpace(void *ida_mem,
  268. long int *lenrwBBDP, long int *leniwBBDP);
  269. SUNDIALS_EXPORT int IDABBDPrecGetNumGfnEvals(void *ida_mem, long int *ngevalsBBDP);
  270. /*
  271. * =================================================================
  272. * PART II - backward problems
  273. * =================================================================
  274. */
  275. /*
  276. * -----------------------------------------------------------------
  277. * Types: IDALocalFnB and IDACommFnB
  278. * -----------------------------------------------------------------
  279. * Local approximation function and inter-process communication
  280. * function for the BBD preconditioner on the backward phase.
  281. * -----------------------------------------------------------------
  282. */
  283. typedef int (*IDABBDLocalFnB)(long int NlocalB, realtype tt,
  284. N_Vector yy, N_Vector yp,
  285. N_Vector yyB, N_Vector ypB, N_Vector gvalB,
  286. void *user_dataB);
  287. typedef int (*IDABBDCommFnB)(long int NlocalB, realtype tt,
  288. N_Vector yy, N_Vector yp,
  289. N_Vector yyB, N_Vector ypB,
  290. void *user_dataB);
  291. /*
  292. * -----------------------------------------------------------------
  293. * Functions: IDABBDPrecInitB, IDABBDPrecReInit
  294. * -----------------------------------------------------------------
  295. * Interface functions for the IDABBDPRE preconditioner to be used
  296. * on the backward phase.
  297. * The 'which' argument is the int returned by IDACreateB.
  298. * -----------------------------------------------------------------
  299. */
  300. SUNDIALS_EXPORT int IDABBDPrecInitB(void *ida_mem, int which, long int NlocalB,
  301. long int mudqB, long int mldqB,
  302. long int mukeepB, long int mlkeepB,
  303. realtype dq_rel_yyB,
  304. IDABBDLocalFnB GresB, IDABBDCommFnB GcommB);
  305. SUNDIALS_EXPORT int IDABBDPrecReInitB(void *ida_mem, int which,
  306. long int mudqB, long int mldqB, realtype dq_rel_yyB);
  307. #ifdef __cplusplus
  308. }
  309. #endif
  310. #endif