Actual source code: gcr.c
1: #include <petsc/private/kspimpl.h>
3: typedef struct {
4: PetscInt restart;
5: PetscInt n_restarts;
6: PetscScalar *val;
7: Vec *VV, *SS;
8: Vec R;
10: KSPFlexibleModifyPCFn *modifypc; /* function to modify the preconditioner*/
11: PetscCtxDestroyFn *modifypc_destroy; /* function to destroy the user context for the modifypc function */
13: void *modifypc_ctx; /* user defined data for the modifypc function */
14: } KSP_GCR;
16: static PetscErrorCode KSPSolve_GCR_cycle(KSP ksp)
17: {
18: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
19: PetscScalar r_dot_v;
20: Mat A, B;
21: PC pc;
22: Vec s, v, r;
23: /*
24: The residual norm will not be computed when ksp->its > ksp->chknorm hence need to initialize norm_r with some dummy value
25: */
26: PetscReal norm_r = 0.0, nrm;
27: PetscInt k, i, restart;
28: Vec x;
30: PetscFunctionBegin;
31: restart = ctx->restart;
32: PetscCall(KSPGetPC(ksp, &pc));
33: PetscCall(KSPGetOperators(ksp, &A, &B));
35: x = ksp->vec_sol;
36: r = ctx->R;
38: for (k = 0; k < restart; k++) {
39: v = ctx->VV[k];
40: s = ctx->SS[k];
41: if (ctx->modifypc) PetscCall((*ctx->modifypc)(ksp, ksp->its, k, ksp->rnorm, ctx->modifypc_ctx));
43: PetscCall(KSP_PCApply(ksp, r, s)); /* s = B^{-1} r */
44: PetscCall(KSP_MatMult(ksp, A, s, v)); /* v = A s */
46: PetscCall(VecMDot(v, k, ctx->VV, ctx->val));
47: for (i = 0; i < k; i++) ctx->val[i] = -ctx->val[i];
48: PetscCall(VecMAXPY(v, k, ctx->val, ctx->VV)); /* v = v - sum_{i=0}^{k-1} alpha_i v_i */
49: PetscCall(VecMAXPY(s, k, ctx->val, ctx->SS)); /* s = s - sum_{i=0}^{k-1} alpha_i s_i */
51: PetscCall(VecDotNorm2(r, v, &r_dot_v, &nrm));
52: nrm = PetscSqrtReal(nrm);
53: r_dot_v = r_dot_v / nrm;
54: PetscCall(VecScale(v, 1.0 / nrm));
55: PetscCall(VecScale(s, 1.0 / nrm));
56: PetscCall(VecAXPY(x, r_dot_v, s));
57: PetscCall(VecAXPY(r, -r_dot_v, v));
58: if (ksp->its > ksp->chknorm && ksp->normtype != KSP_NORM_NONE) {
59: PetscCall(VecNorm(r, NORM_2, &norm_r));
60: KSPCheckNorm(ksp, norm_r);
61: }
62: /* update the local counter and the global counter */
63: ksp->its++;
64: ksp->rnorm = norm_r;
66: PetscCall(KSPLogResidualHistory(ksp, norm_r));
67: PetscCall(KSPMonitor(ksp, ksp->its, norm_r));
69: if (ksp->its - 1 > ksp->chknorm) {
70: PetscCall((*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP));
71: if (ksp->reason) break;
72: }
74: if (ksp->its >= ksp->max_it) {
75: ksp->reason = KSP_CONVERGED_ITS;
76: break;
77: }
78: }
79: ctx->n_restarts++;
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
83: static PetscErrorCode KSPSolve_GCR(KSP ksp)
84: {
85: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
86: Mat A, B;
87: Vec r, b, x;
88: PetscReal norm_r = 0.0;
90: PetscFunctionBegin;
91: PetscCall(KSPGetOperators(ksp, &A, &B));
92: x = ksp->vec_sol;
93: b = ksp->vec_rhs;
94: r = ctx->R;
96: /* compute initial residual */
97: PetscCall(KSP_MatMult(ksp, A, x, r));
98: PetscCall(VecAYPX(r, -1.0, b)); /* r = b - A x */
99: if (ksp->normtype != KSP_NORM_NONE) {
100: PetscCall(VecNorm(r, NORM_2, &norm_r));
101: KSPCheckNorm(ksp, norm_r);
102: }
103: ksp->its = 0;
104: ksp->rnorm0 = norm_r;
106: PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm0));
107: PetscCall(KSPMonitor(ksp, ksp->its, ksp->rnorm0));
108: PetscCall((*ksp->converged)(ksp, ksp->its, ksp->rnorm0, &ksp->reason, ksp->cnvP));
109: if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
111: do {
112: PetscCall(KSPSolve_GCR_cycle(ksp));
113: if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS); /* catch case when convergence occurs inside the cycle */
114: } while (ksp->its < ksp->max_it);
116: if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
117: PetscFunctionReturn(PETSC_SUCCESS);
118: }
120: static PetscErrorCode KSPView_GCR(KSP ksp, PetscViewer viewer)
121: {
122: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
123: PetscBool isascii;
125: PetscFunctionBegin;
126: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
127: if (isascii) {
128: PetscCall(PetscViewerASCIIPrintf(viewer, " restart = %" PetscInt_FMT " \n", ctx->restart));
129: PetscCall(PetscViewerASCIIPrintf(viewer, " restarts performed = %" PetscInt_FMT " \n", ctx->n_restarts));
130: }
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: static PetscErrorCode KSPSetUp_GCR(KSP ksp)
135: {
136: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
137: Mat A;
138: PetscBool diagonalscale;
140: PetscFunctionBegin;
141: PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
142: PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
144: PetscCall(KSPGetOperators(ksp, &A, NULL));
145: PetscCall(MatCreateVecs(A, &ctx->R, NULL));
146: PetscCall(VecDuplicateVecs(ctx->R, ctx->restart, &ctx->VV));
147: PetscCall(VecDuplicateVecs(ctx->R, ctx->restart, &ctx->SS));
149: PetscCall(PetscMalloc1(ctx->restart, &ctx->val));
150: PetscFunctionReturn(PETSC_SUCCESS);
151: }
153: static PetscErrorCode KSPReset_GCR(KSP ksp)
154: {
155: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
157: PetscFunctionBegin;
158: PetscCall(VecDestroy(&ctx->R));
159: PetscCall(VecDestroyVecs(ctx->restart, &ctx->VV));
160: PetscCall(VecDestroyVecs(ctx->restart, &ctx->SS));
161: if (ctx->modifypc_destroy) PetscCall((*ctx->modifypc_destroy)(&ctx->modifypc_ctx));
162: PetscCall(PetscFree(ctx->val));
163: PetscFunctionReturn(PETSC_SUCCESS);
164: }
166: static PetscErrorCode KSPDestroy_GCR(KSP ksp)
167: {
168: PetscFunctionBegin;
169: PetscCall(KSPReset_GCR(ksp));
170: PetscCall(KSPDestroyDefault(ksp));
171: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRSetRestart_C", NULL));
172: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRGetRestart_C", NULL));
173: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRSetModifyPC_C", NULL));
174: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFlexibleSetModifyPC_C", NULL));
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: static PetscErrorCode KSPSetFromOptions_GCR(KSP ksp, PetscOptionItems PetscOptionsObject)
179: {
180: KSP_GCR *ctx = (KSP_GCR *)ksp->data;
181: PetscInt restart;
182: PetscBool flg;
184: PetscFunctionBegin;
185: PetscOptionsHeadBegin(PetscOptionsObject, "KSP GCR options");
186: PetscCall(PetscOptionsInt("-ksp_gcr_restart", "Number of Krylov search directions", "KSPGCRSetRestart", ctx->restart, &restart, &flg));
187: if (flg) PetscCall(KSPGCRSetRestart(ksp, restart));
188: PetscOptionsHeadEnd();
189: PetscFunctionReturn(PETSC_SUCCESS);
190: }
192: static PetscErrorCode KSPGCRSetModifyPC_GCR(KSP ksp, KSPFlexibleModifyPCFn *function, void *ctx, PetscCtxDestroyFn *destroy)
193: {
194: KSP_GCR *gcr = (KSP_GCR *)ksp->data;
196: PetscFunctionBegin;
198: gcr->modifypc = function;
199: gcr->modifypc_destroy = destroy;
200: gcr->modifypc_ctx = ctx;
201: PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: /*@C
205: KSPGCRSetModifyPC - Sets the routine used by `KSPGCR` to modify the preconditioner for each iteration
207: Logically Collective
209: Input Parameters:
210: + ksp - iterative context obtained from `KSPCreate()`
211: . function - user defined function to modify the preconditioner, see `KSPFlexibleModifyPCFn`
212: . ctx - user provided context for the modify preconditioner function
213: - destroy - the function to use to destroy the user provided application context.
215: Level: intermediate
217: .seealso: [](ch_ksp), `KSP`, `KSPGCR`, `KSPFlexibleModifyPCFn`, `KSPFGMRESModifyPCFn`, [](sec_flexibleksp)
218: @*/
219: PetscErrorCode KSPGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *function, void *ctx, PetscCtxDestroyFn *destroy)
220: {
221: PetscFunctionBegin;
222: PetscUseMethod(ksp, "KSPGCRSetModifyPC_C", (KSP, KSPFlexibleModifyPCFn *, void *, PetscCtxDestroyFn *), (ksp, function, ctx, destroy));
223: PetscFunctionReturn(PETSC_SUCCESS);
224: }
226: static PetscErrorCode KSPGCRSetRestart_GCR(KSP ksp, PetscInt restart)
227: {
228: KSP_GCR *ctx;
230: PetscFunctionBegin;
231: ctx = (KSP_GCR *)ksp->data;
232: ctx->restart = restart;
233: PetscFunctionReturn(PETSC_SUCCESS);
234: }
236: static PetscErrorCode KSPGCRGetRestart_GCR(KSP ksp, PetscInt *restart)
237: {
238: KSP_GCR *ctx;
240: PetscFunctionBegin;
241: ctx = (KSP_GCR *)ksp->data;
242: *restart = ctx->restart;
243: PetscFunctionReturn(PETSC_SUCCESS);
244: }
246: /*@
247: KSPGCRSetRestart - Sets number of iterations at which `KSPGCR` restarts.
249: Not Collective
251: Input Parameters:
252: + ksp - the Krylov space context
253: - restart - integer restart value
255: Options Database Key:
256: . -ksp_gcr_restart <restart> - the number of stored vectors to orthogonalize against
258: Level: intermediate
260: Note:
261: The default value is 30.
263: Developer Note:
264: The API could be made uniform for all `KSP` methods that have a restart.
266: .seealso: [](ch_ksp), `KSPGCR`, `KSPSetTolerances()`, `KSPGCRGetRestart()`, `KSPGMRESSetRestart()`
267: @*/
268: PetscErrorCode KSPGCRSetRestart(KSP ksp, PetscInt restart)
269: {
270: PetscFunctionBegin;
271: PetscTryMethod(ksp, "KSPGCRSetRestart_C", (KSP, PetscInt), (ksp, restart));
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: /*@
276: KSPGCRGetRestart - Gets number of iterations at which `KSPGCR` restarts.
278: Not Collective
280: Input Parameter:
281: . ksp - the Krylov space context
283: Output Parameter:
284: . restart - integer restart value
286: Level: intermediate
288: .seealso: [](ch_ksp), `KSPGCR`, `KSPSetTolerances()`, `KSPGCRSetRestart()`, `KSPGMRESGetRestart()`
289: @*/
290: PetscErrorCode KSPGCRGetRestart(KSP ksp, PetscInt *restart)
291: {
292: PetscFunctionBegin;
293: PetscTryMethod(ksp, "KSPGCRGetRestart_C", (KSP, PetscInt *), (ksp, restart));
294: PetscFunctionReturn(PETSC_SUCCESS);
295: }
297: static PetscErrorCode KSPBuildSolution_GCR(KSP ksp, Vec v, Vec *V)
298: {
299: Vec x;
301: PetscFunctionBegin;
302: x = ksp->vec_sol;
303: if (v) {
304: PetscCall(VecCopy(x, v));
305: if (V) *V = v;
306: } else if (V) {
307: *V = ksp->vec_sol;
308: }
309: PetscFunctionReturn(PETSC_SUCCESS);
310: }
312: static PetscErrorCode KSPBuildResidual_GCR(KSP ksp, Vec t, Vec v, Vec *V)
313: {
314: KSP_GCR *ctx;
316: PetscFunctionBegin;
317: ctx = (KSP_GCR *)ksp->data;
318: if (v) {
319: PetscCall(VecCopy(ctx->R, v));
320: if (V) *V = v;
321: } else if (V) {
322: *V = ctx->R;
323: }
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
327: /*MC
328: KSPGCR - Implements the preconditioned flexible Generalized Conjugate Residual method {cite}`eisenstat1983variational`. [](sec_flexibleksp),
330: Options Database Key:
331: . -ksp_gcr_restart <restart> - the number of stored vectors to orthogonalize against
333: Level: beginner
335: Notes:
336: The GCR Krylov method supports non-symmetric matrices and permits the use of a preconditioner
337: which may vary from one iteration to the next.
339: Users can define a method to vary the preconditioner between iterates via `KSPFlexibleSetModifyPC()` or `KSPGCRSetModifyPC()`.
341: Restarts are solves with x0 not equal to zero. When a restart occurs, the initial starting
342: solution is given by the current estimate for x which was obtained by the last restart
343: iterations of the GCR algorithm.
345: Unlike `KSPGMRES` and `KSPFGMRES`, when using GCR, the solution and residual vector can be directly accessed at any iterate,
346: with zero computational cost, via a call to `KSPBuildSolution()` and `KSPBuildResidual()` respectively.
348: This implementation of GCR will only apply the stopping condition test whenever ksp->its > ksp->chknorm,
349: where ksp->chknorm is specified via the command line argument `-ksp_check_norm_iteration` or via
350: the function `KSPSetCheckNormIteration()`. Hence the residual norm reported by the monitor and stored
351: in the residual history will be listed as 0.0 before this iteration. It is actually not 0.0; just not calculated.
353: The method implemented requires the storage of 2 x restart + 1 vectors, twice as much as `KSPGMRES`.
355: Support only for right preconditioning.
357: Contributed by:
358: Dave May
360: .seealso: [](ch_ksp), [](sec_flexibleksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPGCRSetRestart()`, `KSPGCRGetRestart()`,
361: `KSPGCRSetRestart()`, `KSPFlexibleSetModifyPC()`, `KSPGCRSetModifyPC()`, `KSPGMRES`, `KSPFGMRES`
362: M*/
363: PETSC_EXTERN PetscErrorCode KSPCreate_GCR(KSP ksp)
364: {
365: KSP_GCR *ctx;
367: PetscFunctionBegin;
368: PetscCall(PetscNew(&ctx));
370: ctx->restart = 30;
371: ctx->n_restarts = 0;
372: ksp->data = (void *)ctx;
374: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));
375: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 3));
377: ksp->ops->setup = KSPSetUp_GCR;
378: ksp->ops->solve = KSPSolve_GCR;
379: ksp->ops->reset = KSPReset_GCR;
380: ksp->ops->destroy = KSPDestroy_GCR;
381: ksp->ops->view = KSPView_GCR;
382: ksp->ops->setfromoptions = KSPSetFromOptions_GCR;
383: ksp->ops->buildsolution = KSPBuildSolution_GCR;
384: ksp->ops->buildresidual = KSPBuildResidual_GCR;
386: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRSetRestart_C", KSPGCRSetRestart_GCR));
387: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRGetRestart_C", KSPGCRGetRestart_GCR));
388: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRSetModifyPC_C", KSPGCRSetModifyPC_GCR));
389: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFlexibleSetModifyPC_C", KSPGCRSetModifyPC_GCR));
390: PetscFunctionReturn(PETSC_SUCCESS);
391: }