sundials_spgmr.h 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302
  1. /*
  2. * -----------------------------------------------------------------
  3. * $Revision: 4378 $
  4. * $Date: 2015-02-19 10:55:14 -0800 (Thu, 19 Feb 2015) $
  5. * -----------------------------------------------------------------
  6. * Programmer(s): Scott D. Cohen, Alan C. Hindmarsh and
  7. * Radu Serban @ 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 implementation of SPGMR Krylov
  20. * iterative linear solver. The SPGMR algorithm is based on the
  21. * Scaled Preconditioned GMRES (Generalized Minimal Residual)
  22. * method.
  23. *
  24. * The SPGMR algorithm solves a linear system A x = b.
  25. * Preconditioning is allowed on the left, right, or both.
  26. * Scaling is allowed on both sides, and restarts are also allowed.
  27. * We denote the preconditioner and scaling matrices as follows:
  28. * P1 = left preconditioner
  29. * P2 = right preconditioner
  30. * S1 = diagonal matrix of scale factors for P1-inverse b
  31. * S2 = diagonal matrix of scale factors for P2 x
  32. * The matrices A, P1, and P2 are not required explicitly; only
  33. * routines that provide A, P1-inverse, and P2-inverse as
  34. * operators are required.
  35. *
  36. * In this notation, SPGMR applies the underlying GMRES method to
  37. * the equivalent transformed system
  38. * Abar xbar = bbar , where
  39. * Abar = S1 (P1-inverse) A (P2-inverse) (S2-inverse) ,
  40. * bbar = S1 (P1-inverse) b , and xbar = S2 P2 x .
  41. *
  42. * The scaling matrices must be chosen so that vectors S1
  43. * P1-inverse b and S2 P2 x have dimensionless components.
  44. * If preconditioning is done on the left only (P2 = I), by a
  45. * matrix P, then S2 must be a scaling for x, while S1 is a
  46. * scaling for P-inverse b, and so may also be taken as a scaling
  47. * for x. Similarly, if preconditioning is done on the right only
  48. * (P1 = I, P2 = P), then S1 must be a scaling for b, while S2 is
  49. * a scaling for P x, and may also be taken as a scaling for b.
  50. *
  51. * The stopping test for the SPGMR iterations is on the L2 norm of
  52. * the scaled preconditioned residual:
  53. * || bbar - Abar xbar ||_2 < delta
  54. * with an input test constant delta.
  55. *
  56. * The usage of this SPGMR solver involves supplying two routines
  57. * and making three calls. The user-supplied routines are
  58. * atimes (A_data, x, y) to compute y = A x, given x,
  59. * and
  60. * psolve (P_data, y, x, lr)
  61. * to solve P1 x = y or P2 x = y for x, given y.
  62. * The three user calls are:
  63. * mem = SpgmrMalloc(lmax, vec_tmpl);
  64. * to initialize memory,
  65. * flag = SpgmrSolve(mem,A_data,x,b,...,
  66. * P_data,s1,s2,atimes,psolve,...);
  67. * to solve the system, and
  68. * SpgmrFree(mem);
  69. * to free the memory created by SpgmrMalloc.
  70. * Complete details for specifying atimes and psolve and for the
  71. * usage calls are given below and in sundials_iterative.h.
  72. * -----------------------------------------------------------------
  73. */
  74. #ifndef _SPGMR_H
  75. #define _SPGMR_H
  76. #include <sundials/sundials_iterative.h>
  77. #ifdef __cplusplus /* wrapper to enable C++ usage */
  78. extern "C" {
  79. #endif
  80. /*
  81. * -----------------------------------------------------------------
  82. * Types: SpgmrMemRec, SpgmrMem
  83. * -----------------------------------------------------------------
  84. * SpgmrMem is a pointer to an SpgmrMemRec which contains
  85. * the memory needed by SpgmrSolve. The SpgmrMalloc routine
  86. * returns a pointer of type SpgmrMem which should then be passed
  87. * in subsequent calls to SpgmrSolve. The SpgmrFree routine frees
  88. * the memory allocated by SpgmrMalloc.
  89. *
  90. * l_max is the maximum Krylov dimension that SpgmrSolve will be
  91. * permitted to use.
  92. *
  93. * V is the array of Krylov basis vectors v_1, ..., v_(l_max+1),
  94. * stored in V[0], ..., V[l_max], where l_max is the second
  95. * parameter to SpgmrMalloc. Each v_i is a vector of type
  96. * N_Vector.
  97. *
  98. * Hes is the (l_max+1) x l_max Hessenberg matrix. It is stored
  99. * row-wise so that the (i,j)th element is given by Hes[i][j].
  100. *
  101. * givens is a length 2*l_max array which represents the
  102. * Givens rotation matrices that arise in the algorithm. The
  103. * Givens rotation matrices F_0, F_1, ..., F_j, where F_i is
  104. *
  105. * 1
  106. * 1
  107. * c_i -s_i <--- row i
  108. * s_i c_i
  109. * 1
  110. * 1
  111. *
  112. * are represented in the givens vector as
  113. * givens[0]=c_0, givens[1]=s_0, givens[2]=c_1, givens[3]=s_1,
  114. * ..., givens[2j]=c_j, givens[2j+1]=s_j.
  115. *
  116. * xcor is a vector (type N_Vector) which holds the scaled,
  117. * preconditioned correction to the initial guess.
  118. *
  119. * yg is a length (l_max+1) array of realtype used to hold "short"
  120. * vectors (e.g. y and g).
  121. *
  122. * vtemp is a vector (type N_Vector) used as temporary vector
  123. * storage during calculations.
  124. * -----------------------------------------------------------------
  125. */
  126. typedef struct _SpgmrMemRec {
  127. int l_max;
  128. N_Vector *V;
  129. realtype **Hes;
  130. realtype *givens;
  131. N_Vector xcor;
  132. realtype *yg;
  133. N_Vector vtemp;
  134. } SpgmrMemRec, *SpgmrMem;
  135. /*
  136. * -----------------------------------------------------------------
  137. * Function : SpgmrMalloc
  138. * -----------------------------------------------------------------
  139. * SpgmrMalloc allocates the memory used by SpgmrSolve. It
  140. * returns a pointer of type SpgmrMem which the user of the
  141. * SPGMR package should pass to SpgmrSolve. The parameter l_max
  142. * is the maximum Krylov dimension that SpgmrSolve will be
  143. * permitted to use. The parameter vec_tmpl is a pointer to an
  144. * N_Vector used as a template to create new vectors by duplication.
  145. * This routine returns NULL if there is a memory request failure.
  146. * -----------------------------------------------------------------
  147. */
  148. SUNDIALS_EXPORT SpgmrMem SpgmrMalloc(int l_max, N_Vector vec_tmpl);
  149. /*
  150. * -----------------------------------------------------------------
  151. * Function : SpgmrSolve
  152. * -----------------------------------------------------------------
  153. * SpgmrSolve solves the linear system Ax = b using the SPGMR
  154. * method. The return values are given by the symbolic constants
  155. * below. The first SpgmrSolve parameter is a pointer to memory
  156. * allocated by a prior call to SpgmrMalloc.
  157. *
  158. * mem is the pointer returned by SpgmrMalloc to the structure
  159. * containing the memory needed by SpgmrSolve.
  160. *
  161. * A_data is a pointer to information about the coefficient
  162. * matrix A. This pointer is passed to the user-supplied function
  163. * atimes.
  164. *
  165. * x is the initial guess x_0 upon entry and the solution
  166. * N_Vector upon exit with return value SPGMR_SUCCESS or
  167. * SPGMR_RES_REDUCED. For all other return values, the output x
  168. * is undefined.
  169. *
  170. * b is the right hand side N_Vector. It is undisturbed by this
  171. * function.
  172. *
  173. * pretype is the type of preconditioning to be used. Its
  174. * legal values are enumerated in sundials_iterative.h. These
  175. * values are PREC_NONE=0, PREC_LEFT=1, PREC_RIGHT=2, and
  176. * PREC_BOTH=3.
  177. *
  178. * gstype is the type of Gram-Schmidt orthogonalization to be
  179. * used. Its legal values are enumerated in sundials_iterative.h.
  180. * These values are MODIFIED_GS=0 and CLASSICAL_GS=1.
  181. *
  182. * delta is the tolerance on the L2 norm of the scaled,
  183. * preconditioned residual. On return with value SPGMR_SUCCESS,
  184. * this residual satisfies || s1 P1_inv (b - Ax) ||_2 <= delta.
  185. *
  186. * max_restarts is the maximum number of times the algorithm is
  187. * allowed to restart.
  188. *
  189. * P_data is a pointer to preconditioner information. This
  190. * pointer is passed to the user-supplied function psolve.
  191. *
  192. * s1 is an N_Vector of positive scale factors for P1-inv b, where
  193. * P1 is the left preconditioner. (Not tested for positivity.)
  194. * Pass NULL if no scaling on P1-inv b is required.
  195. *
  196. * s2 is an N_Vector of positive scale factors for P2 x, where
  197. * P2 is the right preconditioner. (Not tested for positivity.)
  198. * Pass NULL if no scaling on P2 x is required.
  199. *
  200. * atimes is the user-supplied function which performs the
  201. * operation of multiplying A by a given vector. Its description
  202. * is given in sundials_iterative.h.
  203. *
  204. * psolve is the user-supplied function which solves a
  205. * preconditioner system Pz = r, where P is P1 or P2. Its full
  206. * description is given in sundials_iterative.h. The psolve function
  207. * will not be called if pretype is NONE; in that case, the user
  208. * should pass NULL for psolve.
  209. *
  210. * res_norm is a pointer to the L2 norm of the scaled,
  211. * preconditioned residual. On return with value SPGMR_SUCCESS or
  212. * SPGMR_RES_REDUCED, (*res_norm) contains the value
  213. * || s1 P1_inv (b - Ax) ||_2 for the computed solution x.
  214. * For all other return values, (*res_norm) is undefined. The
  215. * caller is responsible for allocating the memory (*res_norm)
  216. * to be filled in by SpgmrSolve.
  217. *
  218. * nli is a pointer to the number of linear iterations done in
  219. * the execution of SpgmrSolve. The caller is responsible for
  220. * allocating the memory (*nli) to be filled in by SpgmrSolve.
  221. *
  222. * nps is a pointer to the number of calls made to psolve during
  223. * the execution of SpgmrSolve. The caller is responsible for
  224. * allocating the memory (*nps) to be filled in by SpgmrSolve.
  225. *
  226. * Note: Repeated calls can be made to SpgmrSolve with varying
  227. * input arguments. If, however, the problem size N or the
  228. * maximum Krylov dimension l_max changes, then a call to
  229. * SpgmrMalloc must be made to obtain new memory for SpgmrSolve
  230. * to use.
  231. * -----------------------------------------------------------------
  232. */
  233. SUNDIALS_EXPORT int SpgmrSolve(SpgmrMem mem, void *A_data, N_Vector x, N_Vector b,
  234. int pretype, int gstype, realtype delta,
  235. int max_restarts, void *P_data, N_Vector s1,
  236. N_Vector s2, ATimesFn atimes, PSolveFn psolve,
  237. realtype *res_norm, int *nli, int *nps);
  238. /* Return values for SpgmrSolve */
  239. #define SPGMR_SUCCESS 0 /* Converged */
  240. #define SPGMR_RES_REDUCED 1 /* Did not converge, but reduced
  241. norm of residual */
  242. #define SPGMR_CONV_FAIL 2 /* Failed to converge */
  243. #define SPGMR_QRFACT_FAIL 3 /* QRfact found singular matrix */
  244. #define SPGMR_PSOLVE_FAIL_REC 4 /* psolve failed recoverably */
  245. #define SPGMR_ATIMES_FAIL_REC 5 /* atimes failed recoverably */
  246. #define SPGMR_PSET_FAIL_REC 6 /* pset faild recoverably */
  247. #define SPGMR_MEM_NULL -1 /* mem argument is NULL */
  248. #define SPGMR_ATIMES_FAIL_UNREC -2 /* atimes returned failure flag */
  249. #define SPGMR_PSOLVE_FAIL_UNREC -3 /* psolve failed unrecoverably */
  250. #define SPGMR_GS_FAIL -4 /* Gram-Schmidt routine faiuled */
  251. #define SPGMR_QRSOL_FAIL -5 /* QRsol found singular R */
  252. #define SPGMR_PSET_FAIL_UNREC -6 /* pset failed unrecoverably */
  253. /*
  254. * -----------------------------------------------------------------
  255. * Function : SpgmrFree
  256. * -----------------------------------------------------------------
  257. * SpgmrMalloc frees the memory allocated by SpgmrMalloc. It is
  258. * illegal to use the pointer mem after a call to SpgmrFree.
  259. * -----------------------------------------------------------------
  260. */
  261. SUNDIALS_EXPORT void SpgmrFree(SpgmrMem mem);
  262. /*
  263. * -----------------------------------------------------------------
  264. * Macro: SPGMR_VTEMP
  265. * -----------------------------------------------------------------
  266. * This macro provides access to the work vector vtemp in the
  267. * memory block of the SPGMR module. The argument mem is the
  268. * memory pointer returned by SpgmrMalloc, of type SpgmrMem,
  269. * and the macro value is of type N_Vector.
  270. * On a return from SpgmrSolve with *nli = 0, this vector
  271. * contains the scaled preconditioned initial residual,
  272. * s1 * P1_inverse * (b - A x_0).
  273. * -----------------------------------------------------------------
  274. */
  275. #define SPGMR_VTEMP(mem) (mem->vtemp)
  276. #ifdef __cplusplus
  277. }
  278. #endif
  279. #endif