Actual source code: gmres.c
2: /*
3: This file implements GMRES (a Generalized Minimal Residual) method.
4: Reference: Saad and Schultz, 1986.
6: Some comments on left vs. right preconditioning, and restarts.
7: Left and right preconditioning.
8: If right preconditioning is chosen, then the problem being solved
9: by gmres is actually
10: My = AB^-1 y = f
11: so the initial residual is
12: r = f - Mx
13: Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
14: residual is
15: r = f - A x
16: The final solution is then
17: x = B^-1 y
19: If left preconditioning is chosen, then the problem being solved is
20: My = B^-1 A x = B^-1 f,
21: and the initial residual is
22: r = B^-1(f - Ax)
24: Restarts: Restarts are basically solves with x0 not equal to zero.
25: Note that we can eliminate an extra application of B^-1 between
26: restarts as long as we don't require that the solution at the end
27: of an unsuccessful gmres iteration always be the solution x.
28: */
30: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
31: #define GMRES_DELTA_DIRECTIONS 10
32: #define GMRES_DEFAULT_MAXK 30
33: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP, PetscInt, PetscBool, PetscReal *);
34: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);
36: PetscErrorCode KSPSetUp_GMRES(KSP ksp)
37: {
38: PetscInt hh, hes, rs, cc;
39: PetscInt max_k, k;
40: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
42: PetscFunctionBegin;
43: max_k = gmres->max_k; /* restart size */
44: hh = (max_k + 2) * (max_k + 1);
45: hes = (max_k + 1) * (max_k + 1);
46: rs = (max_k + 2);
47: cc = (max_k + 1);
49: PetscCall(PetscCalloc5(hh, &gmres->hh_origin, hes, &gmres->hes_origin, rs, &gmres->rs_origin, cc, &gmres->cc_origin, cc, &gmres->ss_origin));
51: if (ksp->calc_sings) {
52: /* Allocate workspace to hold Hessenberg matrix needed by lapack */
53: PetscCall(PetscMalloc1((max_k + 3) * (max_k + 9), &gmres->Rsvd));
54: PetscCall(PetscMalloc1(6 * (max_k + 2), &gmres->Dsvd));
55: }
57: /* Allocate array to hold pointers to user vectors. Note that we need
58: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
59: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k + gmres->nextra_vecs;
61: PetscCall(PetscMalloc1(gmres->vecs_allocated, &gmres->vecs));
62: PetscCall(PetscMalloc1(VEC_OFFSET + 2 + max_k, &gmres->user_work));
63: PetscCall(PetscMalloc1(VEC_OFFSET + 2 + max_k, &gmres->mwork_alloc));
65: if (gmres->q_preallocate) {
66: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
68: PetscCall(KSPCreateVecs(ksp, gmres->vv_allocated, &gmres->user_work[0], 0, NULL));
70: gmres->mwork_alloc[0] = gmres->vv_allocated;
71: gmres->nwork_alloc = 1;
72: for (k = 0; k < gmres->vv_allocated; k++) gmres->vecs[k] = gmres->user_work[0][k];
73: } else {
74: gmres->vv_allocated = 5;
76: PetscCall(KSPCreateVecs(ksp, 5, &gmres->user_work[0], 0, NULL));
78: gmres->mwork_alloc[0] = 5;
79: gmres->nwork_alloc = 1;
80: for (k = 0; k < gmres->vv_allocated; k++) gmres->vecs[k] = gmres->user_work[0][k];
81: }
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }
85: /*
86: Run gmres, possibly with restart. Return residual history if requested.
87: input parameters:
89: . gmres - structure containing parameters and work areas
91: output parameters:
92: . nres - residuals (from preconditioned system) at each step.
93: If restarting, consider passing nres+it. If null,
94: ignored
95: . itcount - number of iterations used. nres[0] to nres[itcount]
96: are defined. If null, ignored.
98: Notes:
99: On entry, the value in vector VEC_VV(0) should be the initial residual
100: (this allows shortcuts where the initial preconditioned residual is 0).
101: */
102: PetscErrorCode KSPGMRESCycle(PetscInt *itcount, KSP ksp)
103: {
104: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
105: PetscReal res, hapbnd, tt;
106: PetscInt it = 0, max_k = gmres->max_k;
107: PetscBool hapend = PETSC_FALSE;
109: PetscFunctionBegin;
110: if (itcount) *itcount = 0;
111: PetscCall(VecNormalize(VEC_VV(0), &res));
112: KSPCheckNorm(ksp, res);
114: /* the constant .1 is arbitrary, just some measure at how incorrect the residuals are */
115: if ((ksp->rnorm > 0.0) && (PetscAbsReal(res - ksp->rnorm) > gmres->breakdowntol * gmres->rnorm0)) {
116: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_CONV_FAILED, "Residual norm computed by GMRES recursion formula %g is far from the computed residual norm %g at restart, residual norm at start of cycle %g",
117: (double)ksp->rnorm, (double)res, (double)gmres->rnorm0);
118: PetscCall(PetscInfo(ksp, "Residual norm computed by GMRES recursion formula %g is far from the computed residual norm %g at restart, residual norm at start of cycle %g", (double)ksp->rnorm, (double)res, (double)gmres->rnorm0));
119: ksp->reason = KSP_DIVERGED_BREAKDOWN;
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
122: *GRS(0) = gmres->rnorm0 = res;
124: /* check for the convergence */
125: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
126: ksp->rnorm = res;
127: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
128: gmres->it = (it - 1);
129: PetscCall(KSPLogResidualHistory(ksp, res));
130: PetscCall(KSPLogErrorHistory(ksp));
131: PetscCall(KSPMonitor(ksp, ksp->its, res));
132: if (!res) {
133: ksp->reason = KSP_CONVERGED_ATOL;
134: PetscCall(PetscInfo(ksp, "Converged due to zero residual norm on entry\n"));
135: PetscFunctionReturn(PETSC_SUCCESS);
136: }
138: PetscCall((*ksp->converged)(ksp, ksp->its, res, &ksp->reason, ksp->cnvP));
139: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
140: if (it) {
141: PetscCall(KSPLogResidualHistory(ksp, res));
142: PetscCall(KSPLogErrorHistory(ksp));
143: PetscCall(KSPMonitor(ksp, ksp->its, res));
144: }
145: gmres->it = (it - 1);
146: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) PetscCall(KSPGMRESGetNewVectors(ksp, it + 1));
147: PetscCall(KSP_PCApplyBAorAB(ksp, VEC_VV(it), VEC_VV(1 + it), VEC_TEMP_MATOP));
149: /* update hessenberg matrix and do Gram-Schmidt */
150: PetscCall((*gmres->orthog)(ksp, it));
151: if (ksp->reason) break;
153: /* vv(i+1) . vv(i+1) */
154: PetscCall(VecNormalize(VEC_VV(it + 1), &tt));
155: KSPCheckNorm(ksp, tt);
157: /* save the magnitude */
158: *HH(it + 1, it) = tt;
159: *HES(it + 1, it) = tt;
161: /* check for the happy breakdown */
162: hapbnd = PetscAbsScalar(tt / *GRS(it));
163: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
164: if (tt < hapbnd) {
165: PetscCall(PetscInfo(ksp, "Detected happy breakdown, current hapbnd = %14.12e tt = %14.12e\n", (double)hapbnd, (double)tt));
166: hapend = PETSC_TRUE;
167: }
168: PetscCall(KSPGMRESUpdateHessenberg(ksp, it, hapend, &res));
170: it++;
171: gmres->it = (it - 1); /* For converged */
172: ksp->its++;
173: ksp->rnorm = res;
174: if (ksp->reason) break;
176: PetscCall((*ksp->converged)(ksp, ksp->its, res, &ksp->reason, ksp->cnvP));
178: /* Catch error in happy breakdown and signal convergence and break from loop */
179: if (hapend) {
180: if (ksp->normtype == KSP_NORM_NONE) { /* convergence test was skipped in this case */
181: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
182: } else if (!ksp->reason) {
183: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "You reached the happy break down, but convergence was not indicated. Residual norm = %g", (double)res);
184: ksp->reason = KSP_DIVERGED_BREAKDOWN;
185: break;
186: }
187: }
188: }
190: /* Monitor if we know that we will not return for a restart */
191: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
192: PetscCall(KSPLogResidualHistory(ksp, res));
193: PetscCall(KSPLogErrorHistory(ksp));
194: PetscCall(KSPMonitor(ksp, ksp->its, res));
195: }
197: if (itcount) *itcount = it;
199: /*
200: Down here we have to solve for the "best" coefficients of the Krylov
201: columns, add the solution values together, and possibly unwind the
202: preconditioning from the solution
203: */
204: /* Form the solution (or the solution so far) */
205: PetscCall(KSPGMRESBuildSoln(GRS(0), ksp->vec_sol, ksp->vec_sol, ksp, it - 1));
206: PetscFunctionReturn(PETSC_SUCCESS);
207: }
209: PetscErrorCode KSPSolve_GMRES(KSP ksp)
210: {
211: PetscInt its, itcount, i;
212: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
213: PetscBool guess_zero = ksp->guess_zero;
214: PetscInt N = gmres->max_k + 1;
216: PetscFunctionBegin;
217: PetscCheck(!ksp->calc_sings || gmres->Rsvd, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ORDER, "Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
219: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
220: ksp->its = 0;
221: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
223: itcount = 0;
224: gmres->fullcycle = 0;
225: ksp->rnorm = -1.0; /* special marker for KSPGMRESCycle() */
226: while (!ksp->reason || (ksp->rnorm == -1 && ksp->reason == KSP_DIVERGED_PC_FAILED)) {
227: PetscCall(KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs));
228: PetscCall(KSPGMRESCycle(&its, ksp));
229: /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
230: if the cycle is complete for the computation of the Ritz pairs */
231: if (its == gmres->max_k) {
232: gmres->fullcycle++;
233: if (ksp->calc_ritz) {
234: if (!gmres->hes_ritz) {
235: PetscCall(PetscMalloc1(N * N, &gmres->hes_ritz));
236: PetscCall(VecDuplicateVecs(VEC_VV(0), N, &gmres->vecb));
237: }
238: PetscCall(PetscArraycpy(gmres->hes_ritz, gmres->hes_origin, N * N));
239: for (i = 0; i < gmres->max_k + 1; i++) PetscCall(VecCopy(VEC_VV(i), gmres->vecb[i]));
240: }
241: }
242: itcount += its;
243: if (itcount >= ksp->max_it) {
244: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
245: break;
246: }
247: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
248: }
249: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: PetscErrorCode KSPReset_GMRES(KSP ksp)
254: {
255: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
256: PetscInt i;
258: PetscFunctionBegin;
259: /* Free the Hessenberg matrices */
260: PetscCall(PetscFree5(gmres->hh_origin, gmres->hes_origin, gmres->rs_origin, gmres->cc_origin, gmres->ss_origin));
261: PetscCall(PetscFree(gmres->hes_ritz));
263: /* free work vectors */
264: PetscCall(PetscFree(gmres->vecs));
265: for (i = 0; i < gmres->nwork_alloc; i++) PetscCall(VecDestroyVecs(gmres->mwork_alloc[i], &gmres->user_work[i]));
266: gmres->nwork_alloc = 0;
267: if (gmres->vecb) PetscCall(VecDestroyVecs(gmres->max_k + 1, &gmres->vecb));
269: PetscCall(PetscFree(gmres->user_work));
270: PetscCall(PetscFree(gmres->mwork_alloc));
271: PetscCall(PetscFree(gmres->nrs));
272: PetscCall(VecDestroy(&gmres->sol_temp));
273: PetscCall(PetscFree(gmres->Rsvd));
274: PetscCall(PetscFree(gmres->Dsvd));
275: PetscCall(PetscFree(gmres->orthogwork));
277: gmres->vv_allocated = 0;
278: gmres->vecs_allocated = 0;
279: gmres->sol_temp = NULL;
280: PetscFunctionReturn(PETSC_SUCCESS);
281: }
283: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
284: {
285: PetscFunctionBegin;
286: PetscCall(KSPReset_GMRES(ksp));
287: PetscCall(PetscFree(ksp->data));
288: /* clear composed functions */
289: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", NULL));
290: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", NULL));
291: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", NULL));
292: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", NULL));
293: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", NULL));
294: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", NULL));
295: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetBreakdownTolerance_C", NULL));
296: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", NULL));
297: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", NULL));
298: PetscFunctionReturn(PETSC_SUCCESS);
299: }
300: /*
301: KSPGMRESBuildSoln - create the solution from the starting vector and the
302: current iterates.
304: Input parameters:
305: nrs - work area of size it + 1.
306: vs - index of initial guess
307: vdest - index of result. Note that vs may == vdest (replace
308: guess with the solution).
310: This is an internal routine that knows about the GMRES internals.
311: */
312: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs, Vec vs, Vec vdest, KSP ksp, PetscInt it)
313: {
314: PetscScalar tt;
315: PetscInt ii, k, j;
316: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
318: PetscFunctionBegin;
319: /* Solve for solution vector that minimizes the residual */
321: /* If it is < 0, no gmres steps have been performed */
322: if (it < 0) {
323: PetscCall(VecCopy(vs, vdest)); /* VecCopy() is smart, exists immediately if vguess == vdest */
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
326: if (*HH(it, it) != 0.0) {
327: nrs[it] = *GRS(it) / *HH(it, it);
328: } else {
329: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "You reached the break down in GMRES; HH(it,it) = 0");
330: ksp->reason = KSP_DIVERGED_BREAKDOWN;
332: PetscCall(PetscInfo(ksp, "Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %" PetscInt_FMT " GRS(it) = %g\n", it, (double)PetscAbsScalar(*GRS(it))));
333: PetscFunctionReturn(PETSC_SUCCESS);
334: }
335: for (ii = 1; ii <= it; ii++) {
336: k = it - ii;
337: tt = *GRS(k);
338: for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
339: if (*HH(k, k) == 0.0) {
340: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %" PetscInt_FMT, k);
341: ksp->reason = KSP_DIVERGED_BREAKDOWN;
342: PetscCall(PetscInfo(ksp, "Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %" PetscInt_FMT "\n", k));
343: PetscFunctionReturn(PETSC_SUCCESS);
344: }
345: nrs[k] = tt / *HH(k, k);
346: }
348: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
349: PetscCall(VecSet(VEC_TEMP, 0.0));
350: PetscCall(VecMAXPY(VEC_TEMP, it + 1, nrs, &VEC_VV(0)));
352: PetscCall(KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP));
353: /* add solution to previous solution */
354: if (vdest != vs) PetscCall(VecCopy(vs, vdest));
355: PetscCall(VecAXPY(vdest, 1.0, VEC_TEMP));
356: PetscFunctionReturn(PETSC_SUCCESS);
357: }
358: /*
359: Do the scalar work for the orthogonalization. Return new residual norm.
360: */
361: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
362: {
363: PetscScalar *hh, *cc, *ss, tt;
364: PetscInt j;
365: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
367: PetscFunctionBegin;
368: hh = HH(0, it);
369: cc = CC(0);
370: ss = SS(0);
372: /* Apply all the previously computed plane rotations to the new column
373: of the Hessenberg matrix */
374: for (j = 1; j <= it; j++) {
375: tt = *hh;
376: *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
377: hh++;
378: *hh = *cc++ * *hh - (*ss++ * tt);
379: }
381: /*
382: compute the new plane rotation, and apply it to:
383: 1) the right-hand-side of the Hessenberg system
384: 2) the new column of the Hessenberg matrix
385: thus obtaining the updated value of the residual
386: */
387: if (!hapend) {
388: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
389: if (tt == 0.0) {
390: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "tt == 0.0");
391: ksp->reason = KSP_DIVERGED_NULL;
392: PetscFunctionReturn(PETSC_SUCCESS);
393: }
394: *cc = *hh / tt;
395: *ss = *(hh + 1) / tt;
396: *GRS(it + 1) = -(*ss * *GRS(it));
397: *GRS(it) = PetscConj(*cc) * *GRS(it);
398: *hh = PetscConj(*cc) * *hh + *ss * *(hh + 1);
399: *res = PetscAbsScalar(*GRS(it + 1));
400: } else {
401: /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
402: another rotation matrix (so RH doesn't change). The new residual is
403: always the new sine term times the residual from last time (GRS(it)),
404: but now the new sine rotation would be zero...so the residual should
405: be zero...so we will multiply "zero" by the last residual. This might
406: not be exactly what we want to do here -could just return "zero". */
408: *res = 0.0;
409: }
410: PetscFunctionReturn(PETSC_SUCCESS);
411: }
412: /*
413: This routine allocates more work vectors, starting from VEC_VV(it).
414: */
415: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp, PetscInt it)
416: {
417: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
418: PetscInt nwork = gmres->nwork_alloc, k, nalloc;
420: PetscFunctionBegin;
421: nalloc = PetscMin(ksp->max_it, gmres->delta_allocate);
422: /* Adjust the number to allocate to make sure that we don't exceed the
423: number of available slots */
424: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
425: if (!nalloc) PetscFunctionReturn(PETSC_SUCCESS);
427: gmres->vv_allocated += nalloc;
429: PetscCall(KSPCreateVecs(ksp, nalloc, &gmres->user_work[nwork], 0, NULL));
431: gmres->mwork_alloc[nwork] = nalloc;
432: for (k = 0; k < nalloc; k++) gmres->vecs[it + VEC_OFFSET + k] = gmres->user_work[nwork][k];
433: gmres->nwork_alloc++;
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }
437: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp, Vec ptr, Vec *result)
438: {
439: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
441: PetscFunctionBegin;
442: if (!ptr) {
443: if (!gmres->sol_temp) PetscCall(VecDuplicate(ksp->vec_sol, &gmres->sol_temp));
444: ptr = gmres->sol_temp;
445: }
446: if (!gmres->nrs) {
447: /* allocate the work area */
448: PetscCall(PetscMalloc1(gmres->max_k, &gmres->nrs));
449: }
451: PetscCall(KSPGMRESBuildSoln(gmres->nrs, ksp->vec_sol, ptr, ksp, gmres->it));
452: if (result) *result = ptr;
453: PetscFunctionReturn(PETSC_SUCCESS);
454: }
456: PetscErrorCode KSPView_GMRES(KSP ksp, PetscViewer viewer)
457: {
458: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
459: const char *cstr;
460: PetscBool iascii, isstring;
462: PetscFunctionBegin;
463: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
464: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
465: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
466: switch (gmres->cgstype) {
467: case (KSP_GMRES_CGS_REFINE_NEVER):
468: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
469: break;
470: case (KSP_GMRES_CGS_REFINE_ALWAYS):
471: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
472: break;
473: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
474: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
475: break;
476: default:
477: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Unknown orthogonalization");
478: }
479: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
480: cstr = "Modified Gram-Schmidt Orthogonalization";
481: } else {
482: cstr = "unknown orthogonalization";
483: }
484: if (iascii) {
485: PetscCall(PetscViewerASCIIPrintf(viewer, " restart=%" PetscInt_FMT ", using %s\n", gmres->max_k, cstr));
486: PetscCall(PetscViewerASCIIPrintf(viewer, " happy breakdown tolerance %g\n", (double)gmres->haptol));
487: } else if (isstring) {
488: PetscCall(PetscViewerStringSPrintf(viewer, "%s restart %" PetscInt_FMT, cstr, gmres->max_k));
489: }
490: PetscFunctionReturn(PETSC_SUCCESS);
491: }
493: /*@C
494: KSPGMRESMonitorKrylov - Calls `VecView()` for each new direction in the `KSPGMRES` accumulated Krylov space.
496: Collective
498: Input Parameters:
499: + ksp - the `KSP` context
500: . its - iteration number
501: . fgnorm - 2-norm of residual (or gradient)
502: - dummy - an collection of viewers created with `KSPViewerCreate()`
504: Options Database Key:
505: . -ksp_gmres_krylov_monitor <bool> - Plot the Krylov directions
507: Level: intermediate
509: Note:
510: A new `PETSCVIEWERDRAW` is created for each Krylov vector so they can all be simultaneously viewed
512: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `VecView()`, `KSPViewersCreate()`, `KSPViewersDestroy()`
513: @*/
514: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp, PetscInt its, PetscReal fgnorm, void *dummy)
515: {
516: PetscViewers viewers = (PetscViewers)dummy;
517: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
518: Vec x;
519: PetscViewer viewer;
520: PetscBool flg;
522: PetscFunctionBegin;
523: PetscCall(PetscViewersGetViewer(viewers, gmres->it + 1, &viewer));
524: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &flg));
525: if (!flg) {
526: PetscCall(PetscViewerSetType(viewer, PETSCVIEWERDRAW));
527: PetscCall(PetscViewerDrawSetInfo(viewer, NULL, "Krylov GMRES Monitor", PETSC_DECIDE, PETSC_DECIDE, 300, 300));
528: }
529: x = VEC_VV(gmres->it + 1);
530: PetscCall(VecView(x, viewer));
531: PetscFunctionReturn(PETSC_SUCCESS);
532: }
534: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
535: {
536: PetscInt restart;
537: PetscReal haptol, breakdowntol;
538: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
539: PetscBool flg;
541: PetscFunctionBegin;
542: PetscOptionsHeadBegin(PetscOptionsObject, "KSP GMRES Options");
543: PetscCall(PetscOptionsInt("-ksp_gmres_restart", "Number of Krylov search directions", "KSPGMRESSetRestart", gmres->max_k, &restart, &flg));
544: if (flg) PetscCall(KSPGMRESSetRestart(ksp, restart));
545: PetscCall(PetscOptionsReal("-ksp_gmres_haptol", "Tolerance for exact convergence (happy ending)", "KSPGMRESSetHapTol", gmres->haptol, &haptol, &flg));
546: if (flg) PetscCall(KSPGMRESSetHapTol(ksp, haptol));
547: PetscCall(PetscOptionsReal("-ksp_gmres_breakdown_tolerance", "Divergence breakdown tolerance during GMRES restart", "KSPGMRESSetBreakdownTolerance", gmres->breakdowntol, &breakdowntol, &flg));
548: if (flg) PetscCall(KSPGMRESSetBreakdownTolerance(ksp, breakdowntol));
549: flg = PETSC_FALSE;
550: PetscCall(PetscOptionsBool("-ksp_gmres_preallocate", "Preallocate Krylov vectors", "KSPGMRESSetPreAllocateVectors", flg, &flg, NULL));
551: if (flg) PetscCall(KSPGMRESSetPreAllocateVectors(ksp));
552: PetscCall(PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt", "Classical (unmodified) Gram-Schmidt (fast)", "KSPGMRESSetOrthogonalization", &flg));
553: if (flg) PetscCall(KSPGMRESSetOrthogonalization(ksp, KSPGMRESClassicalGramSchmidtOrthogonalization));
554: PetscCall(PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt", "Modified Gram-Schmidt (slow,more stable)", "KSPGMRESSetOrthogonalization", &flg));
555: if (flg) PetscCall(KSPGMRESSetOrthogonalization(ksp, KSPGMRESModifiedGramSchmidtOrthogonalization));
556: PetscCall(PetscOptionsEnum("-ksp_gmres_cgs_refinement_type", "Type of iterative refinement for classical (unmodified) Gram-Schmidt", "KSPGMRESSetCGSRefinementType", KSPGMRESCGSRefinementTypes, (PetscEnum)gmres->cgstype, (PetscEnum *)&gmres->cgstype, &flg));
557: flg = PETSC_FALSE;
558: PetscCall(PetscOptionsBool("-ksp_gmres_krylov_monitor", "Plot the Krylov directions", "KSPMonitorSet", flg, &flg, NULL));
559: if (flg) {
560: PetscViewers viewers;
561: PetscCall(PetscViewersCreate(PetscObjectComm((PetscObject)ksp), &viewers));
562: PetscCall(KSPMonitorSet(ksp, KSPGMRESMonitorKrylov, viewers, (PetscErrorCode(*)(void **))PetscViewersDestroy));
563: }
564: PetscOptionsHeadEnd();
565: PetscFunctionReturn(PETSC_SUCCESS);
566: }
568: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp, PetscReal tol)
569: {
570: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
572: PetscFunctionBegin;
573: PetscCheck(tol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Tolerance must be non-negative");
574: gmres->haptol = tol;
575: PetscFunctionReturn(PETSC_SUCCESS);
576: }
578: PetscErrorCode KSPGMRESSetBreakdownTolerance_GMRES(KSP ksp, PetscReal tol)
579: {
580: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
582: PetscFunctionBegin;
583: if (tol == PETSC_DEFAULT) {
584: gmres->breakdowntol = 0.1;
585: PetscFunctionReturn(PETSC_SUCCESS);
586: }
587: PetscCheck(tol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Breakdown tolerance must be non-negative");
588: gmres->breakdowntol = tol;
589: PetscFunctionReturn(PETSC_SUCCESS);
590: }
592: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp, PetscInt *max_k)
593: {
594: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
596: PetscFunctionBegin;
597: *max_k = gmres->max_k;
598: PetscFunctionReturn(PETSC_SUCCESS);
599: }
601: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp, PetscInt max_k)
602: {
603: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
605: PetscFunctionBegin;
606: PetscCheck(max_k >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Restart must be positive");
607: if (!ksp->setupstage) {
608: gmres->max_k = max_k;
609: } else if (gmres->max_k != max_k) {
610: gmres->max_k = max_k;
611: ksp->setupstage = KSP_SETUP_NEW;
612: /* free the data structures, then create them again */
613: PetscCall(KSPReset_GMRES(ksp));
614: }
615: PetscFunctionReturn(PETSC_SUCCESS);
616: }
618: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp, FCN fcn)
619: {
620: PetscFunctionBegin;
621: ((KSP_GMRES *)ksp->data)->orthog = fcn;
622: PetscFunctionReturn(PETSC_SUCCESS);
623: }
625: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp, FCN *fcn)
626: {
627: PetscFunctionBegin;
628: *fcn = ((KSP_GMRES *)ksp->data)->orthog;
629: PetscFunctionReturn(PETSC_SUCCESS);
630: }
632: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
633: {
634: KSP_GMRES *gmres;
636: PetscFunctionBegin;
637: gmres = (KSP_GMRES *)ksp->data;
638: gmres->q_preallocate = 1;
639: PetscFunctionReturn(PETSC_SUCCESS);
640: }
642: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp, KSPGMRESCGSRefinementType type)
643: {
644: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
646: PetscFunctionBegin;
647: gmres->cgstype = type;
648: PetscFunctionReturn(PETSC_SUCCESS);
649: }
651: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp, KSPGMRESCGSRefinementType *type)
652: {
653: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
655: PetscFunctionBegin;
656: *type = gmres->cgstype;
657: PetscFunctionReturn(PETSC_SUCCESS);
658: }
660: /*@
661: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
662: in the classical Gram Schmidt orthogonalization.
664: Logically Collective
666: Input Parameters:
667: + ksp - the Krylov space context
668: - type - the type of refinement
669: .vb
670: KSP_GMRES_CGS_REFINE_NEVER
671: KSP_GMRES_CGS_REFINE_IFNEEDED
672: KSP_GMRES_CGS_REFINE_ALWAYS
673: .ve
675: Options Database Key:
676: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - refinement type
678: Level: intermediate
680: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESGetCGSRefinementType()`,
681: `KSPGMRESGetOrthogonalization()`
682: @*/
683: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp, KSPGMRESCGSRefinementType type)
684: {
685: PetscFunctionBegin;
688: PetscTryMethod(ksp, "KSPGMRESSetCGSRefinementType_C", (KSP, KSPGMRESCGSRefinementType), (ksp, type));
689: PetscFunctionReturn(PETSC_SUCCESS);
690: }
692: /*@
693: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
694: in the classical Gram Schmidt orthogonalization.
696: Not Collective
698: Input Parameter:
699: . ksp - the Krylov space context
701: Output Parameter:
702: . type - the type of refinement
704: Options Database Key:
705: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - type of refinement
707: Level: intermediate
709: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetCGSRefinementType()`,
710: `KSPGMRESGetOrthogonalization()`
711: @*/
712: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp, KSPGMRESCGSRefinementType *type)
713: {
714: PetscFunctionBegin;
716: PetscUseMethod(ksp, "KSPGMRESGetCGSRefinementType_C", (KSP, KSPGMRESCGSRefinementType *), (ksp, type));
717: PetscFunctionReturn(PETSC_SUCCESS);
718: }
720: /*@
721: KSPGMRESSetRestart - Sets number of iterations at which `KSPGMRES`, `KSPFGMRES` and `KSPLGMRES` restarts.
723: Logically Collective
725: Input Parameters:
726: + ksp - the Krylov space context
727: - restart - integer restart value
729: Options Database Key:
730: . -ksp_gmres_restart <positive integer> - integer restart value
732: Level: intermediate
734: Note:
735: The default value is 30.
737: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESGetRestart()`
738: @*/
739: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
740: {
741: PetscFunctionBegin;
744: PetscTryMethod(ksp, "KSPGMRESSetRestart_C", (KSP, PetscInt), (ksp, restart));
745: PetscFunctionReturn(PETSC_SUCCESS);
746: }
748: /*@
749: KSPGMRESGetRestart - Gets number of iterations at which `KSPGMRES`, `KSPFGMRES` and `KSPLGMRES` restarts.
751: Not Collective
753: Input Parameter:
754: . ksp - the Krylov space context
756: Output Parameter:
757: . restart - integer restart value
759: Level: intermediate
761: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetRestart()`
762: @*/
763: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
764: {
765: PetscFunctionBegin;
766: PetscUseMethod(ksp, "KSPGMRESGetRestart_C", (KSP, PetscInt *), (ksp, restart));
767: PetscFunctionReturn(PETSC_SUCCESS);
768: }
770: /*@
771: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in `KSPGMRES`, `KSPFGMRES` and `KSPLGMRES`
773: Logically Collective
775: Input Parameters:
776: + ksp - the Krylov space context
777: - tol - the tolerance
779: Options Database Key:
780: . -ksp_gmres_haptol <positive real value> - set tolerance for determining happy breakdown
782: Level: intermediate
784: Note:
785: Happy breakdown is the rare case in `KSPGMRES` where an 'exact' solution is obtained after
786: a certain number of iterations. If you attempt more iterations after this point unstable
787: things can happen hence very occasionally you may need to set this value to detect this condition
789: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPSetTolerances()`
790: @*/
791: PetscErrorCode KSPGMRESSetHapTol(KSP ksp, PetscReal tol)
792: {
793: PetscFunctionBegin;
795: PetscTryMethod((ksp), "KSPGMRESSetHapTol_C", (KSP, PetscReal), ((ksp), (tol)));
796: PetscFunctionReturn(PETSC_SUCCESS);
797: }
799: /*@
800: KSPGMRESSetBreakdownTolerance - Sets tolerance for determining divergence breakdown in `KSPGMRES`.
802: Logically Collective
804: Input Parameters:
805: + ksp - the Krylov space context
806: - tol - the tolerance
808: Options Database Key:
809: . -ksp_gmres_breakdown_tolerance <positive real value> - set tolerance for determining divergence breakdown
811: Level: intermediate
813: Note:
814: Divergence breakdown occurs when GMRES residual increases significantly during restart
816: .seealso: [](chapter_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetHapTol()`
817: @*/
818: PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP ksp, PetscReal tol)
819: {
820: PetscFunctionBegin;
822: PetscTryMethod((ksp), "KSPGMRESSetBreakdownTolerance_C", (KSP, PetscReal), (ksp, tol));
823: PetscFunctionReturn(PETSC_SUCCESS);
824: }
826: /*MC
827: KSPGMRES - Implements the Generalized Minimal Residual method [1] with restart
829: Options Database Keys:
830: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
831: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
832: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
833: vectors are allocated as needed)
834: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
835: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
836: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
837: stability of the classical Gram-Schmidt orthogonalization.
838: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
840: Level: beginner
842: Note:
843: Left and right preconditioning are supported, but not symmetric preconditioning.
845: Reference:
846: . [1] - YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
847: SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.
849: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`,
850: `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
851: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
852: `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`
853: M*/
855: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
856: {
857: KSP_GMRES *gmres;
859: PetscFunctionBegin;
860: PetscCall(PetscNew(&gmres));
861: ksp->data = (void *)gmres;
863: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 4));
864: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 3));
865: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_SYMMETRIC, 2));
866: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));
867: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
869: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
870: ksp->ops->setup = KSPSetUp_GMRES;
871: ksp->ops->solve = KSPSolve_GMRES;
872: ksp->ops->reset = KSPReset_GMRES;
873: ksp->ops->destroy = KSPDestroy_GMRES;
874: ksp->ops->view = KSPView_GMRES;
875: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
876: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
877: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
878: ksp->ops->computeritz = KSPComputeRitz_GMRES;
879: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES));
880: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES));
881: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", KSPGMRESGetOrthogonalization_GMRES));
882: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES));
883: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_GMRES));
884: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES));
885: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetBreakdownTolerance_C", KSPGMRESSetBreakdownTolerance_GMRES));
886: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES));
887: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", KSPGMRESGetCGSRefinementType_GMRES));
889: gmres->haptol = 1.0e-30;
890: gmres->breakdowntol = 0.1;
891: gmres->q_preallocate = 0;
892: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
893: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
894: gmres->nrs = NULL;
895: gmres->sol_temp = NULL;
896: gmres->max_k = GMRES_DEFAULT_MAXK;
897: gmres->Rsvd = NULL;
898: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
899: gmres->orthogwork = NULL;
900: PetscFunctionReturn(PETSC_SUCCESS);
901: }