1 | /* |
---|
2 | * ----------------------------------------------------------------- |
---|
3 | * Programmer(s): Daniel Reynolds @ SMU |
---|
4 | * Based on code sundials_spgmr.h by: Scott D. Cohen, |
---|
5 | * Alan C. Hindmarsh and Radu Serban @ LLNL |
---|
6 | * ----------------------------------------------------------------- |
---|
7 | * LLNS/SMU Copyright Start |
---|
8 | * Copyright (c) 2017, Southern Methodist University and |
---|
9 | * Lawrence Livermore National Security |
---|
10 | * |
---|
11 | * This work was performed under the auspices of the U.S. Department |
---|
12 | * of Energy by Southern Methodist University and Lawrence Livermore |
---|
13 | * National Laboratory under Contract DE-AC52-07NA27344. |
---|
14 | * Produced at Southern Methodist University and the Lawrence |
---|
15 | * Livermore National Laboratory. |
---|
16 | * |
---|
17 | * All rights reserved. |
---|
18 | * For details, see the LICENSE file. |
---|
19 | * LLNS/SMU Copyright End |
---|
20 | * ----------------------------------------------------------------- |
---|
21 | * This is the header file for the SPGMR implementation of the |
---|
22 | * SUNLINSOL module. The SPGMR algorithm is based on the |
---|
23 | * Scaled Preconditioned GMRES (Generalized Minimal Residual) |
---|
24 | * method. |
---|
25 | * |
---|
26 | * The SPGMR algorithm solves a linear system A x = b. |
---|
27 | * Preconditioning is allowed on the left, right, or both. |
---|
28 | * Scaling is allowed on both sides, and restarts are also allowed. |
---|
29 | * We denote the preconditioner and scaling matrices as follows: |
---|
30 | * P1 = left preconditioner |
---|
31 | * P2 = right preconditioner |
---|
32 | * S1 = diagonal matrix of scale factors for P1-inverse b |
---|
33 | * S2 = diagonal matrix of scale factors for P2 x |
---|
34 | * The matrices A, P1, and P2 are not required explicitly; only |
---|
35 | * routines that provide A, P1-inverse, and P2-inverse as |
---|
36 | * operators are required. |
---|
37 | * |
---|
38 | * In this notation, SPGMR applies the underlying GMRES method to |
---|
39 | * the equivalent transformed system |
---|
40 | * Abar xbar = bbar , where |
---|
41 | * Abar = S1 (P1-inverse) A (P2-inverse) (S2-inverse) , |
---|
42 | * bbar = S1 (P1-inverse) b , and xbar = S2 P2 x . |
---|
43 | * |
---|
44 | * The scaling matrices must be chosen so that vectors S1 |
---|
45 | * P1-inverse b and S2 P2 x have dimensionless components. |
---|
46 | * If preconditioning is done on the left only (P2 = I), by a |
---|
47 | * matrix P, then S2 must be a scaling for x, while S1 is a |
---|
48 | * scaling for P-inverse b, and so may also be taken as a scaling |
---|
49 | * for x. Similarly, if preconditioning is done on the right only |
---|
50 | * (P1 = I, P2 = P), then S1 must be a scaling for b, while S2 is |
---|
51 | * a scaling for P x, and may also be taken as a scaling for b. |
---|
52 | * |
---|
53 | * The stopping test for the SPGMR iterations is on the L2 norm of |
---|
54 | * the scaled preconditioned residual: |
---|
55 | * || bbar - Abar xbar ||_2 < delta |
---|
56 | * with an input test constant delta. |
---|
57 | * |
---|
58 | * The usage of this SPGMR solver involves supplying three routines |
---|
59 | * and making a variety of calls. The user-supplied routines are |
---|
60 | * atimes (A_data, x, y) to compute y = A x, given x, |
---|
61 | * psolve (P_data, y, x, lr) to solve P1 x = y or P2 x = y for |
---|
62 | * x, given y, |
---|
63 | * psetup (P_data) to perform any 'setup' operations in |
---|
64 | * preparation for calling psolve. |
---|
65 | * The user calls are: |
---|
66 | * SUNLinearSolver LS = SUNSPGMR(y, pretype, maxl); |
---|
67 | * to create the linear solver structure, |
---|
68 | * flag = SUNLinSolSetATimes(LS, A_data, atimes); |
---|
69 | * to set the matrix-vector product setup/apply routines, |
---|
70 | * flag = SUNLinSolSetPreconditioner(LS, P_data, psetup, psolve); |
---|
71 | * to *optionally* set the preconditioner setup/apply routines, |
---|
72 | * flag = SUNLinSolInitialize(LS); |
---|
73 | * to perform internal solver memory allocations, |
---|
74 | * flag = SUNLinSolSetup(LS, NULL); |
---|
75 | * to call the psetup routine (if non-NULL); |
---|
76 | * flag = SUNLinSolSolve(LS, NULL, x, b, w, tol); |
---|
77 | * to solve the linear system to the tolerance 'tol' |
---|
78 | * long int nli = SUNLinSolNumIters(LS); |
---|
79 | * to *optionally* retrieve the number of linear iterations |
---|
80 | * performed by the solver, |
---|
81 | * long int lastflag = SUNLinSolLastFlag(LS); |
---|
82 | * to *optionally* retrieve the last internal solver error flag, |
---|
83 | * flag = SUNLinSolFree(LS); |
---|
84 | * to free the solver memory. |
---|
85 | * Complete details for specifying atimes, psetup and psolve |
---|
86 | * and for the usage calls are given below. |
---|
87 | * |
---|
88 | * ----------------------------------------------------------------- |
---|
89 | * |
---|
90 | * Part I contains declarations specific to the SPGMR implementation |
---|
91 | * of the supplied SUNLINSOL module. |
---|
92 | * |
---|
93 | * Part II contains the prototype for the constructor |
---|
94 | * SUNSPGMR as well as implementation-specific prototypes |
---|
95 | * for various useful solver operations. |
---|
96 | * |
---|
97 | * Notes: |
---|
98 | * |
---|
99 | * - The definition of the generic SUNLinearSolver structure can |
---|
100 | * be found in the header file sundials_linearsolver.h. |
---|
101 | * |
---|
102 | * ----------------------------------------------------------------- |
---|
103 | */ |
---|
104 | |
---|
105 | #ifndef _SUNLINSOL_SPGMR_H |
---|
106 | #define _SUNLINSOL_SPGMR_H |
---|
107 | |
---|
108 | #include <sundials/sundials_linearsolver.h> |
---|
109 | #include <sundials/sundials_matrix.h> |
---|
110 | #include <sundials/sundials_nvector.h> |
---|
111 | #include <sundials/sundials_spgmr.h> |
---|
112 | |
---|
113 | #ifdef __cplusplus /* wrapper to enable C++ usage */ |
---|
114 | extern "C" { |
---|
115 | #endif |
---|
116 | |
---|
117 | /* Default SPGMR solver parameters */ |
---|
118 | #define SUNSPGMR_MAXL_DEFAULT 5 |
---|
119 | #define SUNSPGMR_MAXRS_DEFAULT 0 |
---|
120 | #define SUNSPGMR_GSTYPE_DEFAULT MODIFIED_GS |
---|
121 | |
---|
122 | /* |
---|
123 | * ----------------------------------------------------------------- |
---|
124 | * PART I: SPGMR implementation of SUNLinearSolver |
---|
125 | * |
---|
126 | * The SPGMR implementation of the SUNLinearSolver 'content' |
---|
127 | * structure contains: |
---|
128 | * maxl -- number of GMRES basis vectors to use |
---|
129 | * pretype -- flag for type of preconditioning to employ |
---|
130 | * gstype -- flag for type of Gram-Schmidt orthogonalization |
---|
131 | * max_restarts -- number of GMRES restarts to allow |
---|
132 | * numiters -- number of iterations from most-recent solve |
---|
133 | * resnorm -- final linear residual norm from most-recent solve |
---|
134 | * last_flag -- last error return flag from internal setup/solve |
---|
135 | * ATimes -- function pointer to ATimes routine |
---|
136 | * ATData -- pointer to structure for ATimes |
---|
137 | * Psetup -- function pointer to preconditioner setup routine |
---|
138 | * Psolve -- function pointer to preconditioner solve routine |
---|
139 | * PData -- pointer to structure for Psetup/Psolve |
---|
140 | * s1, s2 -- vector pointers for supplied scaling matrices |
---|
141 | * V -- the array of Krylov basis vectors v_1, ..., v_(maxl+1), |
---|
142 | * stored in V[0], ..., V[maxl]. Each v_i is a vector of |
---|
143 | * type N_Vector. |
---|
144 | * Hes -- the (maxl+1) x maxl Hessenberg matrix. It is stored |
---|
145 | * row-wise so that the (i,j)th element is given by Hes[i][j]. |
---|
146 | * givens -- a length 2*max array which represents the Givens |
---|
147 | * rotation matrices that arise in the algorithm. The Givens |
---|
148 | * rotation matrices F_0, F_1, ..., F_j, where F_i is |
---|
149 | * |
---|
150 | * 1 |
---|
151 | * 1 |
---|
152 | * c_i -s_i <--- row i |
---|
153 | * s_i c_i |
---|
154 | * 1 |
---|
155 | * 1 |
---|
156 | * |
---|
157 | * are represented in the givens vector as |
---|
158 | * givens[0]=c_0, givens[1]=s_0, givens[2]=c_1, givens[3]=s_1, |
---|
159 | * ..., givens[2j]=c_j, givens[2j+1]=s_j. |
---|
160 | * xcor -- a vector (type N_Vector) which holds the scaled, |
---|
161 | * preconditioned correction to the initial guess |
---|
162 | * yg -- a length (maxl+1) array of realtype used to hold "short" |
---|
163 | * vectors (e.g. y and g). |
---|
164 | * vtemp -- a vector used as temporary vector storage |
---|
165 | * ----------------------------------------------------------------- |
---|
166 | */ |
---|
167 | |
---|
168 | struct _SUNLinearSolverContent_SPGMR { |
---|
169 | int maxl; |
---|
170 | int pretype; |
---|
171 | int gstype; |
---|
172 | int max_restarts; |
---|
173 | int numiters; |
---|
174 | realtype resnorm; |
---|
175 | long int last_flag; |
---|
176 | |
---|
177 | ATimesFn ATimes; |
---|
178 | void* ATData; |
---|
179 | PSetupFn Psetup; |
---|
180 | PSolveFn Psolve; |
---|
181 | void* PData; |
---|
182 | |
---|
183 | N_Vector s1; |
---|
184 | N_Vector s2; |
---|
185 | N_Vector *V; |
---|
186 | realtype **Hes; |
---|
187 | realtype *givens; |
---|
188 | N_Vector xcor; |
---|
189 | realtype *yg; |
---|
190 | N_Vector vtemp; |
---|
191 | }; |
---|
192 | |
---|
193 | typedef struct _SUNLinearSolverContent_SPGMR *SUNLinearSolverContent_SPGMR; |
---|
194 | |
---|
195 | /* |
---|
196 | * ----------------------------------------------------------------- |
---|
197 | * PART II: functions exported by sunlinsol_spgmr |
---|
198 | * |
---|
199 | * CONSTRUCTOR: |
---|
200 | * SUNSPGMR creates and allocates memory for a SPGMR solver |
---|
201 | * |
---|
202 | * "SET" ROUTINES: |
---|
203 | * SUNSPGMRSetPrecType updates the type of preconditioning to |
---|
204 | * use. Supported values are PREC_NONE, PREC_LEFT, PREC_RIGHT |
---|
205 | * and PREC_BOTH. |
---|
206 | * SUNSPGMRSetGSType sets the type of Gram-Schmidt |
---|
207 | * orthogonalization to use. Supported values are MODIFIED_GS |
---|
208 | * and CLASSICAL_GS. |
---|
209 | * SUNSPGMRSetMaxRestarts sets the number of GMRES restarts to |
---|
210 | * allow. A negative input will result in the default of 0. |
---|
211 | * |
---|
212 | * ----------------------------------------------------------------- |
---|
213 | */ |
---|
214 | |
---|
215 | SUNDIALS_EXPORT SUNLinearSolver SUNSPGMR(N_Vector y, int pretype, int maxl); |
---|
216 | SUNDIALS_EXPORT int SUNSPGMRSetPrecType(SUNLinearSolver S, int pretype); |
---|
217 | SUNDIALS_EXPORT int SUNSPGMRSetGSType(SUNLinearSolver S, int gstype); |
---|
218 | SUNDIALS_EXPORT int SUNSPGMRSetMaxRestarts(SUNLinearSolver S, int maxrs); |
---|
219 | |
---|
220 | /* |
---|
221 | * ----------------------------------------------------------------- |
---|
222 | * SPGMR implementations of various useful linear solver operations |
---|
223 | * ----------------------------------------------------------------- |
---|
224 | */ |
---|
225 | |
---|
226 | SUNDIALS_EXPORT SUNLinearSolver_Type SUNLinSolGetType_SPGMR(SUNLinearSolver S); |
---|
227 | SUNDIALS_EXPORT int SUNLinSolInitialize_SPGMR(SUNLinearSolver S); |
---|
228 | SUNDIALS_EXPORT int SUNLinSolSetATimes_SPGMR(SUNLinearSolver S, void* A_data, |
---|
229 | ATimesFn ATimes); |
---|
230 | SUNDIALS_EXPORT int SUNLinSolSetPreconditioner_SPGMR(SUNLinearSolver S, |
---|
231 | void* P_data, |
---|
232 | PSetupFn Pset, |
---|
233 | PSolveFn Psol); |
---|
234 | SUNDIALS_EXPORT int SUNLinSolSetScalingVectors_SPGMR(SUNLinearSolver S, |
---|
235 | N_Vector s1, |
---|
236 | N_Vector s2); |
---|
237 | SUNDIALS_EXPORT int SUNLinSolSetup_SPGMR(SUNLinearSolver S, SUNMatrix A); |
---|
238 | SUNDIALS_EXPORT int SUNLinSolSolve_SPGMR(SUNLinearSolver S, SUNMatrix A, |
---|
239 | N_Vector x, N_Vector b, realtype tol); |
---|
240 | SUNDIALS_EXPORT int SUNLinSolNumIters_SPGMR(SUNLinearSolver S); |
---|
241 | SUNDIALS_EXPORT realtype SUNLinSolResNorm_SPGMR(SUNLinearSolver S); |
---|
242 | SUNDIALS_EXPORT N_Vector SUNLinSolResid_SPGMR(SUNLinearSolver S); |
---|
243 | SUNDIALS_EXPORT long int SUNLinSolLastFlag_SPGMR(SUNLinearSolver S); |
---|
244 | SUNDIALS_EXPORT int SUNLinSolSpace_SPGMR(SUNLinearSolver S, |
---|
245 | long int *lenrwLS, |
---|
246 | long int *leniwLS); |
---|
247 | SUNDIALS_EXPORT int SUNLinSolFree_SPGMR(SUNLinearSolver S); |
---|
248 | |
---|
249 | |
---|
250 | #ifdef __cplusplus |
---|
251 | } |
---|
252 | #endif |
---|
253 | |
---|
254 | #endif |
---|