Actual source code: gmres.c
1: /*
2: This file implements GMRES (a Generalized Minimal Residual) method.
3: Reference: Saad and Schultz, 1986.
5: Some comments on left vs. right preconditioning, and restarts.
6: Left and right preconditioning.
7: If right preconditioning is chosen, then the problem being solved
8: by GMRES is actually
9: My = AB^-1 y = f
10: so the initial residual is
11: r = f - M y
12: Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
13: residual is
14: r = f - A x
15: The final solution is then
16: x = B^-1 y
18: If left preconditioning is chosen, then the problem being solved is
19: My = B^-1 A x = B^-1 f,
20: and the initial residual is
21: r = B^-1(f - Ax)
23: Restarts: Restarts are basically solves with x0 not equal to zero.
24: Note that we can eliminate an extra application of B^-1 between
25: restarts as long as we don't require that the solution at the end
26: of an unsuccessful gmres iteration always be the solution x.
27: */
29: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
30: #define GMRES_DELTA_DIRECTIONS 10
31: #define GMRES_DEFAULT_MAXK 30
32: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP, PetscInt, PetscBool, PetscReal *);
33: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);
35: PetscErrorCode KSPSetUp_GMRES(KSP ksp)
36: {
37: PetscInt hh, hes, rs, cc;
38: PetscInt max_k, k;
39: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
41: PetscFunctionBegin;
42: max_k = gmres->max_k; /* restart size */
43: hh = (max_k + 2) * (max_k + 1);
44: hes = (max_k + 1) * (max_k + 1);
45: rs = (max_k + 2);
46: cc = (max_k + 1);
48: PetscCall(PetscCalloc5(hh, &gmres->hh_origin, hes, &gmres->hes_origin, rs, &gmres->rs_origin, cc, &gmres->cc_origin, cc, &gmres->ss_origin));
50: if (ksp->calc_sings) {
51: /* Allocate workspace to hold Hessenberg matrix needed by LAPACK */
52: PetscCall(PetscMalloc1((max_k + 3) * (max_k + 9), &gmres->Rsvd));
53: PetscCall(PetscMalloc1(6 * (max_k + 2), &gmres->Dsvd));
54: }
56: /* Allocate array to hold pointers to user vectors. Note that we need
57: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
58: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k + gmres->nextra_vecs;
59: PetscCall(PetscMalloc1(gmres->vecs_allocated, &gmres->vecs));
60: PetscCall(PetscMalloc1(VEC_OFFSET + 2 + max_k, &gmres->user_work));
61: PetscCall(PetscMalloc1(VEC_OFFSET + 2 + max_k, &gmres->mwork_alloc));
62: if (gmres->q_preallocate || ksp->normtype == KSP_NORM_NONE) gmres->vv_allocated = VEC_OFFSET + 2 + PetscMin(max_k, ksp->max_it);
63: else gmres->vv_allocated = VEC_OFFSET + 2 + PetscMin(PetscMin(5, max_k), ksp->max_it);
64: PetscCall(KSPCreateVecs(ksp, gmres->vv_allocated, &gmres->user_work[0], 0, NULL));
65: gmres->mwork_alloc[0] = gmres->vv_allocated;
66: gmres->nwork_alloc = 1;
67: for (k = 0; k < gmres->vv_allocated; k++) gmres->vecs[k] = gmres->user_work[0][k];
68: PetscFunctionReturn(PETSC_SUCCESS);
69: }
71: /*
72: Run gmres, possibly with restart. Return residual history if requested.
73: input parameters:
75: . gmres - structure containing parameters and work areas
77: output parameters:
78: . nres - residuals (from preconditioned system) at each step.
79: If restarting, consider passing nres+it. If null,
80: ignored
81: . itcount - number of iterations used. nres[0] to nres[itcount]
82: are defined. If null, ignored.
84: Notes:
85: On entry, the value in vector VEC_VV(0) should be the initial residual
86: (this allows shortcuts where the initial preconditioned residual is 0).
87: */
88: static PetscErrorCode KSPGMRESCycle(PetscInt *itcount, KSP ksp)
89: {
90: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
91: PetscReal res, hapbnd, tt;
92: PetscInt it = 0, max_k = gmres->max_k;
93: PetscBool hapend = PETSC_FALSE;
95: PetscFunctionBegin;
96: if (itcount) *itcount = 0;
97: PetscCall(VecNormalize(VEC_VV(0), &res));
98: KSPCheckNorm(ksp, res);
100: /* the constant .1 is arbitrary, just some measure at how incorrect the residuals are */
101: if ((ksp->rnorm > 0.0) && (PetscAbsReal(res - ksp->rnorm) > gmres->breakdowntol * gmres->rnorm0)) {
102: 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",
103: (double)ksp->rnorm, (double)res, (double)gmres->rnorm0);
104: 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\n", (double)ksp->rnorm, (double)res, (double)gmres->rnorm0));
105: ksp->reason = KSP_DIVERGED_BREAKDOWN;
106: PetscFunctionReturn(PETSC_SUCCESS);
107: }
108: *GRS(0) = gmres->rnorm0 = res;
110: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
111: ksp->rnorm = res;
112: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
113: gmres->it = (it - 1);
114: PetscCall(KSPLogResidualHistory(ksp, res));
115: PetscCall(KSPLogErrorHistory(ksp));
116: PetscCall(KSPMonitor(ksp, ksp->its, res));
117: if (!res) {
118: ksp->reason = KSP_CONVERGED_ATOL;
119: PetscCall(PetscInfo(ksp, "Converged due to zero residual norm on entry\n"));
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
123: /* check for the convergence */
124: PetscCall((*ksp->converged)(ksp, ksp->its, res, &ksp->reason, ksp->cnvP));
125: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
126: if (it) {
127: PetscCall(KSPLogResidualHistory(ksp, res));
128: PetscCall(KSPLogErrorHistory(ksp));
129: PetscCall(KSPMonitor(ksp, ksp->its, res));
130: }
131: gmres->it = (it - 1);
132: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) PetscCall(KSPGMRESGetNewVectors(ksp, it + 1));
133: PetscCall(KSP_PCApplyBAorAB(ksp, VEC_VV(it), VEC_VV(1 + it), VEC_TEMP_MATOP));
135: /* update Hessenberg matrix and do Gram-Schmidt */
136: PetscCall((*gmres->orthog)(ksp, it));
137: if (ksp->reason) break;
139: /* vv(i+1) . vv(i+1) */
140: PetscCall(VecNormalize(VEC_VV(it + 1), &tt));
141: KSPCheckNorm(ksp, tt);
143: /* save the magnitude */
144: *HH(it + 1, it) = tt;
145: *HES(it + 1, it) = tt;
147: /* check for the happy breakdown */
148: hapbnd = PetscAbsScalar(tt / *GRS(it));
149: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
150: if (tt < hapbnd) {
151: PetscCall(PetscInfo(ksp, "Detected happy ending, current hapbnd = %14.12e tt = %14.12e\n", (double)hapbnd, (double)tt));
152: hapend = PETSC_TRUE;
153: }
154: PetscCall(KSPGMRESUpdateHessenberg(ksp, it, hapend, &res));
156: it++;
157: gmres->it = (it - 1); /* For converged */
158: ksp->its++;
159: ksp->rnorm = res;
160: if (ksp->reason) break;
162: PetscCall((*ksp->converged)(ksp, ksp->its, res, &ksp->reason, ksp->cnvP));
164: /* Catch error in happy breakdown and signal convergence and break from loop */
165: if (hapend) {
166: if (ksp->normtype == KSP_NORM_NONE) { /* convergence test was skipped in this case */
167: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
168: } else if (!ksp->reason) {
169: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Reached happy break down, but convergence was not indicated. Residual norm = %g", (double)res);
170: ksp->reason = KSP_DIVERGED_BREAKDOWN;
171: break;
172: }
173: }
174: }
176: if (itcount) *itcount = it;
178: /*
179: Down here we have to solve for the "best" coefficients of the Krylov
180: columns, add the solution values together, and possibly unwind the
181: preconditioning from the solution
182: */
183: /* Form the solution (or the solution so far) */
184: PetscCall(KSPGMRESBuildSoln(GRS(0), ksp->vec_sol, ksp->vec_sol, ksp, it - 1));
186: /* Monitor if we know that we will not return for a restart */
187: if (ksp->reason == KSP_CONVERGED_ITERATING && ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
188: if (it && ksp->reason) {
189: PetscCall(KSPLogResidualHistory(ksp, res));
190: PetscCall(KSPLogErrorHistory(ksp));
191: PetscCall(KSPMonitor(ksp, ksp->its, res));
192: }
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: static PetscErrorCode KSPSolve_GMRES(KSP ksp)
197: {
198: PetscInt its, itcount, i;
199: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
200: PetscBool guess_zero = ksp->guess_zero;
201: PetscInt N = gmres->max_k + 1;
203: PetscFunctionBegin;
204: PetscCheck(!ksp->calc_sings || gmres->Rsvd, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ORDER, "Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
206: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
207: ksp->its = 0;
208: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
210: itcount = 0;
211: gmres->fullcycle = 0;
212: ksp->rnorm = -1.0; /* special marker for KSPGMRESCycle() */
213: while (!ksp->reason || (ksp->rnorm == -1 && ksp->reason == KSP_DIVERGED_PC_FAILED)) {
214: PetscCall(KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs));
215: PetscCall(KSPGMRESCycle(&its, ksp));
216: /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
217: if the cycle is complete for the computation of the Ritz pairs */
218: if (its == gmres->max_k) {
219: gmres->fullcycle++;
220: if (ksp->calc_ritz) {
221: if (!gmres->hes_ritz) {
222: PetscCall(PetscMalloc1(N * N, &gmres->hes_ritz));
223: PetscCall(VecDuplicateVecs(VEC_VV(0), N, &gmres->vecb));
224: }
225: PetscCall(PetscArraycpy(gmres->hes_ritz, gmres->hes_origin, N * N));
226: for (i = 0; i < gmres->max_k + 1; i++) PetscCall(VecCopy(VEC_VV(i), gmres->vecb[i]));
227: }
228: }
229: itcount += its;
230: if (itcount >= ksp->max_it) {
231: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
232: break;
233: }
234: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
235: }
236: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
237: PetscFunctionReturn(PETSC_SUCCESS);
238: }
240: PetscErrorCode KSPReset_GMRES(KSP ksp)
241: {
242: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
243: PetscInt i;
245: PetscFunctionBegin;
246: /* Free the Hessenberg matrices */
247: PetscCall(PetscFree5(gmres->hh_origin, gmres->hes_origin, gmres->rs_origin, gmres->cc_origin, gmres->ss_origin));
248: PetscCall(PetscFree(gmres->hes_ritz));
250: /* free work vectors */
251: PetscCall(PetscFree(gmres->vecs));
252: for (i = 0; i < gmres->nwork_alloc; i++) PetscCall(VecDestroyVecs(gmres->mwork_alloc[i], &gmres->user_work[i]));
253: gmres->nwork_alloc = 0;
254: if (gmres->vecb) PetscCall(VecDestroyVecs(gmres->max_k + 1, &gmres->vecb));
256: PetscCall(PetscFree(gmres->user_work));
257: PetscCall(PetscFree(gmres->mwork_alloc));
258: PetscCall(PetscFree(gmres->nrs));
259: PetscCall(VecDestroy(&gmres->sol_temp));
260: PetscCall(PetscFree(gmres->Rsvd));
261: PetscCall(PetscFree(gmres->Dsvd));
262: PetscCall(PetscFree(gmres->orthogwork));
264: gmres->vv_allocated = 0;
265: gmres->vecs_allocated = 0;
266: gmres->sol_temp = NULL;
267: PetscFunctionReturn(PETSC_SUCCESS);
268: }
270: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
271: {
272: PetscFunctionBegin;
273: PetscCall(KSPReset_GMRES(ksp));
274: PetscCall(PetscFree(ksp->data));
275: /* clear composed functions */
276: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", NULL));
277: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", NULL));
278: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", NULL));
279: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", NULL));
280: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", NULL));
281: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", NULL));
282: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetBreakdownTolerance_C", NULL));
283: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", NULL));
284: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", NULL));
285: PetscFunctionReturn(PETSC_SUCCESS);
286: }
287: /*
288: KSPGMRESBuildSoln - create the solution from the starting vector and the
289: current iterates.
291: Input parameters:
292: nrs - work area of size it + 1.
293: vs - index of initial guess
294: vdest - index of result. Note that vs may == vdest (replace
295: guess with the solution).
297: This is an internal routine that knows about the GMRES internals.
298: */
299: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs, Vec vs, Vec vdest, KSP ksp, PetscInt it)
300: {
301: PetscScalar tt;
302: PetscInt ii, k, j;
303: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
305: PetscFunctionBegin;
306: /* Solve for solution vector that minimizes the residual */
308: /* If it is < 0, no gmres steps have been performed */
309: if (it < 0) {
310: PetscCall(VecCopy(vs, vdest)); /* VecCopy() is smart, exists immediately if vguess == vdest */
311: PetscFunctionReturn(PETSC_SUCCESS);
312: }
313: if (*HH(it, it) != 0.0) {
314: nrs[it] = *GRS(it) / *HH(it, it);
315: } else {
316: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "You reached the break down in GMRES; HH(it,it) = 0");
317: ksp->reason = KSP_DIVERGED_BREAKDOWN;
319: 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))));
320: PetscFunctionReturn(PETSC_SUCCESS);
321: }
322: for (ii = 1; ii <= it; ii++) {
323: k = it - ii;
324: tt = *GRS(k);
325: for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
326: if (*HH(k, k) == 0.0) {
327: 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);
328: ksp->reason = KSP_DIVERGED_BREAKDOWN;
329: PetscCall(PetscInfo(ksp, "Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %" PetscInt_FMT "\n", k));
330: PetscFunctionReturn(PETSC_SUCCESS);
331: }
332: nrs[k] = tt / *HH(k, k);
333: }
335: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
336: PetscCall(VecMAXPBY(VEC_TEMP, it + 1, nrs, 0, &VEC_VV(0)));
338: PetscCall(KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP));
339: /* add solution to previous solution */
340: if (vdest != vs) PetscCall(VecCopy(vs, vdest));
341: PetscCall(VecAXPY(vdest, 1.0, VEC_TEMP));
342: PetscFunctionReturn(PETSC_SUCCESS);
343: }
344: /*
345: Do the scalar work for the orthogonalization. Return new residual norm.
346: */
347: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
348: {
349: PetscScalar *hh, *cc, *ss, tt;
350: PetscInt j;
351: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
353: PetscFunctionBegin;
354: hh = HH(0, it);
355: cc = CC(0);
356: ss = SS(0);
358: /* Apply all the previously computed plane rotations to the new column
359: of the Hessenberg matrix */
360: for (j = 1; j <= it; j++) {
361: tt = *hh;
362: *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
363: hh++;
364: *hh = *cc++ * *hh - (*ss++ * tt);
365: }
367: /*
368: compute the new plane rotation, and apply it to:
369: 1) the right-hand side of the Hessenberg system
370: 2) the new column of the Hessenberg matrix
371: thus obtaining the updated value of the residual
372: */
373: if (!hapend) {
374: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
375: if (tt == 0.0) {
376: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "tt == 0.0");
377: ksp->reason = KSP_DIVERGED_NULL;
378: PetscFunctionReturn(PETSC_SUCCESS);
379: }
380: *cc = *hh / tt;
381: *ss = *(hh + 1) / tt;
382: *GRS(it + 1) = -(*ss * *GRS(it));
383: *GRS(it) = PetscConj(*cc) * *GRS(it);
384: *hh = PetscConj(*cc) * *hh + *ss * *(hh + 1);
385: *res = PetscAbsScalar(*GRS(it + 1));
386: } else {
387: /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
388: another rotation matrix (so RH doesn't change). The new residual is
389: always the new sine term times the residual from last time (GRS(it)),
390: but now the new sine rotation would be zero...so the residual should
391: be zero...so we will multiply "zero" by the last residual. This might
392: not be exactly what we want to do here -could just return "zero". */
394: *res = 0.0;
395: }
396: PetscFunctionReturn(PETSC_SUCCESS);
397: }
398: /*
399: This routine allocates more work vectors, starting from VEC_VV(it).
400: */
401: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp, PetscInt it)
402: {
403: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
404: PetscInt nwork = gmres->nwork_alloc, k, nalloc;
406: PetscFunctionBegin;
407: nalloc = PetscMin(ksp->max_it, gmres->delta_allocate);
408: /* Adjust the number to allocate to make sure that we don't exceed the
409: number of available slots */
410: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
411: if (!nalloc) PetscFunctionReturn(PETSC_SUCCESS);
413: gmres->vv_allocated += nalloc;
415: PetscCall(KSPCreateVecs(ksp, nalloc, &gmres->user_work[nwork], 0, NULL));
417: gmres->mwork_alloc[nwork] = nalloc;
418: for (k = 0; k < nalloc; k++) gmres->vecs[it + VEC_OFFSET + k] = gmres->user_work[nwork][k];
419: gmres->nwork_alloc++;
420: PetscFunctionReturn(PETSC_SUCCESS);
421: }
423: static PetscErrorCode KSPBuildSolution_GMRES(KSP ksp, Vec ptr, Vec *result)
424: {
425: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
427: PetscFunctionBegin;
428: if (!ptr) {
429: if (!gmres->sol_temp) PetscCall(VecDuplicate(ksp->vec_sol, &gmres->sol_temp));
430: ptr = gmres->sol_temp;
431: }
432: if (!gmres->nrs) {
433: /* allocate the work area */
434: PetscCall(PetscMalloc1(gmres->max_k, &gmres->nrs));
435: }
437: PetscCall(KSPGMRESBuildSoln(gmres->nrs, ksp->vec_sol, ptr, ksp, gmres->it));
438: if (result) *result = ptr;
439: PetscFunctionReturn(PETSC_SUCCESS);
440: }
442: PetscErrorCode KSPView_GMRES(KSP ksp, PetscViewer viewer)
443: {
444: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
445: const char *cstr;
446: PetscBool isascii, isstring;
448: PetscFunctionBegin;
449: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
450: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
451: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
452: switch (gmres->cgstype) {
453: case (KSP_GMRES_CGS_REFINE_NEVER):
454: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
455: break;
456: case (KSP_GMRES_CGS_REFINE_ALWAYS):
457: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
458: break;
459: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
460: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
461: break;
462: default:
463: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Unknown orthogonalization");
464: }
465: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
466: cstr = "Modified Gram-Schmidt Orthogonalization";
467: } else {
468: cstr = "unknown orthogonalization";
469: }
470: if (isascii) {
471: PetscCall(PetscViewerASCIIPrintf(viewer, " restart=%" PetscInt_FMT ", using %s\n", gmres->max_k, cstr));
472: PetscCall(PetscViewerASCIIPrintf(viewer, " happy breakdown tolerance %g\n", (double)gmres->haptol));
473: } else if (isstring) {
474: PetscCall(PetscViewerStringSPrintf(viewer, "%s restart %" PetscInt_FMT, cstr, gmres->max_k));
475: }
476: PetscFunctionReturn(PETSC_SUCCESS);
477: }
479: /*@C
480: KSPGMRESMonitorKrylov - Calls `VecView()` to monitor each new direction in the `KSPGMRES` accumulated Krylov space.
482: Collective
484: Input Parameters:
485: + ksp - the `KSP` context
486: . its - iteration number
487: . fgnorm - 2-norm of residual (or gradient)
488: - dummy - a collection of viewers created with `PetscViewersCreate()`
490: Options Database Key:
491: . -ksp_gmres_krylov_monitor <bool> - Plot the Krylov directions
493: Level: intermediate
495: Note:
496: A new `PETSCVIEWERDRAW` is created for each Krylov vector so they can all be simultaneously viewed
498: .seealso: [](ch_ksp), `KSPGMRES`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `VecView()`, `PetscViewersCreate()`, `PetscViewersDestroy()`
499: @*/
500: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp, PetscInt its, PetscReal fgnorm, void *dummy)
501: {
502: PetscViewers viewers = (PetscViewers)dummy;
503: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
504: Vec x;
505: PetscViewer viewer;
506: PetscBool flg;
508: PetscFunctionBegin;
509: PetscCall(PetscViewersGetViewer(viewers, gmres->it + 1, &viewer));
510: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &flg));
511: if (!flg) {
512: PetscCall(PetscViewerSetType(viewer, PETSCVIEWERDRAW));
513: PetscCall(PetscViewerDrawSetInfo(viewer, NULL, "Krylov GMRES Monitor", PETSC_DECIDE, PETSC_DECIDE, 300, 300));
514: }
515: x = VEC_VV(gmres->it + 1);
516: PetscCall(VecView(x, viewer));
517: PetscFunctionReturn(PETSC_SUCCESS);
518: }
520: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp, PetscOptionItems PetscOptionsObject)
521: {
522: PetscInt restart;
523: PetscReal haptol, breakdowntol;
524: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
525: PetscBool flg;
527: PetscFunctionBegin;
528: PetscOptionsHeadBegin(PetscOptionsObject, "KSP GMRES Options");
529: PetscCall(PetscOptionsInt("-ksp_gmres_restart", "Number of Krylov search directions", "KSPGMRESSetRestart", gmres->max_k, &restart, &flg));
530: if (flg) PetscCall(KSPGMRESSetRestart(ksp, restart));
531: PetscCall(PetscOptionsReal("-ksp_gmres_haptol", "Tolerance for exact convergence (happy ending)", "KSPGMRESSetHapTol", gmres->haptol, &haptol, &flg));
532: if (flg) PetscCall(KSPGMRESSetHapTol(ksp, haptol));
533: PetscCall(PetscOptionsReal("-ksp_gmres_breakdown_tolerance", "Divergence breakdown tolerance during GMRES restart", "KSPGMRESSetBreakdownTolerance", gmres->breakdowntol, &breakdowntol, &flg));
534: if (flg) PetscCall(KSPGMRESSetBreakdownTolerance(ksp, breakdowntol));
535: flg = PETSC_FALSE;
536: PetscCall(PetscOptionsBool("-ksp_gmres_preallocate", "Preallocate Krylov vectors", "KSPGMRESSetPreAllocateVectors", flg, &flg, NULL));
537: if (flg) PetscCall(KSPGMRESSetPreAllocateVectors(ksp));
538: PetscCall(PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt", "Classical (unmodified) Gram-Schmidt (fast)", "KSPGMRESSetOrthogonalization", &flg));
539: if (flg) PetscCall(KSPGMRESSetOrthogonalization(ksp, KSPGMRESClassicalGramSchmidtOrthogonalization));
540: PetscCall(PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt", "Modified Gram-Schmidt (slow,more stable)", "KSPGMRESSetOrthogonalization", &flg));
541: if (flg) PetscCall(KSPGMRESSetOrthogonalization(ksp, KSPGMRESModifiedGramSchmidtOrthogonalization));
542: 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));
543: flg = PETSC_FALSE;
544: PetscCall(PetscOptionsBool("-ksp_gmres_krylov_monitor", "Plot the Krylov directions", "KSPMonitorSet", flg, &flg, NULL));
545: if (flg) {
546: PetscViewers viewers;
548: PetscCall(PetscViewersCreate(PetscObjectComm((PetscObject)ksp), &viewers));
549: PetscCall(KSPMonitorSet(ksp, KSPGMRESMonitorKrylov, viewers, (PetscCtxDestroyFn *)PetscViewersDestroy));
550: }
551: PetscOptionsHeadEnd();
552: PetscFunctionReturn(PETSC_SUCCESS);
553: }
555: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp, PetscReal tol)
556: {
557: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
559: PetscFunctionBegin;
560: PetscCheck(tol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Tolerance must be non-negative");
561: gmres->haptol = tol;
562: PetscFunctionReturn(PETSC_SUCCESS);
563: }
565: static PetscErrorCode KSPGMRESSetBreakdownTolerance_GMRES(KSP ksp, PetscReal tol)
566: {
567: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
569: PetscFunctionBegin;
570: if (tol == (PetscReal)PETSC_DEFAULT) {
571: gmres->breakdowntol = 0.1;
572: PetscFunctionReturn(PETSC_SUCCESS);
573: }
574: PetscCheck(tol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Breakdown tolerance must be non-negative");
575: gmres->breakdowntol = tol;
576: PetscFunctionReturn(PETSC_SUCCESS);
577: }
579: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp, PetscInt *max_k)
580: {
581: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
583: PetscFunctionBegin;
584: *max_k = gmres->max_k;
585: PetscFunctionReturn(PETSC_SUCCESS);
586: }
588: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp, PetscInt max_k)
589: {
590: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
592: PetscFunctionBegin;
593: PetscCheck(max_k >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Restart must be positive");
594: if (!ksp->setupstage) {
595: gmres->max_k = max_k;
596: } else if (gmres->max_k != max_k) {
597: gmres->max_k = max_k;
598: ksp->setupstage = KSP_SETUP_NEW;
599: /* free the data structures, then create them again */
600: PetscCall(KSPReset_GMRES(ksp));
601: }
602: PetscFunctionReturn(PETSC_SUCCESS);
603: }
605: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp, FCN fcn)
606: {
607: PetscFunctionBegin;
608: ((KSP_GMRES *)ksp->data)->orthog = fcn;
609: PetscFunctionReturn(PETSC_SUCCESS);
610: }
612: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp, FCN *fcn)
613: {
614: PetscFunctionBegin;
615: *fcn = ((KSP_GMRES *)ksp->data)->orthog;
616: PetscFunctionReturn(PETSC_SUCCESS);
617: }
619: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
620: {
621: KSP_GMRES *gmres;
623: PetscFunctionBegin;
624: gmres = (KSP_GMRES *)ksp->data;
625: gmres->q_preallocate = 1;
626: PetscFunctionReturn(PETSC_SUCCESS);
627: }
629: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp, KSPGMRESCGSRefinementType type)
630: {
631: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
633: PetscFunctionBegin;
634: gmres->cgstype = type;
635: PetscFunctionReturn(PETSC_SUCCESS);
636: }
638: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp, KSPGMRESCGSRefinementType *type)
639: {
640: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
642: PetscFunctionBegin;
643: *type = gmres->cgstype;
644: PetscFunctionReturn(PETSC_SUCCESS);
645: }
647: /*@
648: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
649: in the classical Gram-Schmidt orthogonalization used by `KSPGMRES` and other PETSc GMRES implementations.
651: Logically Collective
653: Input Parameters:
654: + ksp - the Krylov space solver context
655: - type - the type of refinement
656: .vb
657: KSP_GMRES_CGS_REFINE_NEVER
658: KSP_GMRES_CGS_REFINE_IFNEEDED
659: KSP_GMRES_CGS_REFINE_ALWAYS
660: .ve
662: Options Database Key:
663: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - refinement type
665: Level: intermediate
667: Notes:
668: The default is `KSP_GMRES_CGS_REFINE_NEVER`
670: For a very small set of problems not using refinement, that is `KSP_GMRES_CGS_REFINE_NEVER` may be unstable, thus causing `KSPSolve()`
671: to not converge.
673: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESGetCGSRefinementType()`,
674: `KSPGMRESGetOrthogonalization()`
675: @*/
676: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp, KSPGMRESCGSRefinementType type)
677: {
678: PetscFunctionBegin;
681: PetscTryMethod(ksp, "KSPGMRESSetCGSRefinementType_C", (KSP, KSPGMRESCGSRefinementType), (ksp, type));
682: PetscFunctionReturn(PETSC_SUCCESS);
683: }
685: /*@
686: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
687: in the classical Gram-Schmidt orthogonalization used by `KSPGMRES` and other PETSc GMRES implementations.
689: Not Collective
691: Input Parameter:
692: . ksp - the Krylov space solver context
694: Output Parameter:
695: . type - the type of refinement
697: Level: intermediate
699: .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetCGSRefinementType()`,
700: `KSPGMRESGetOrthogonalization()`
701: @*/
702: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp, KSPGMRESCGSRefinementType *type)
703: {
704: PetscFunctionBegin;
706: PetscUseMethod(ksp, "KSPGMRESGetCGSRefinementType_C", (KSP, KSPGMRESCGSRefinementType *), (ksp, type));
707: PetscFunctionReturn(PETSC_SUCCESS);
708: }
710: /*@
711: KSPGMRESSetRestart - Sets number of iterations at which GMRES (`KSPGMRES`, `KSPFGMRES`, `KSPPGMRES`, `KSPAGMRES`, `KSPDGMRES`, `KSPPIPEFGMRES`,
712: and `KSPLGMRES`) restarts.
714: Logically Collective
716: Input Parameters:
717: + ksp - the Krylov space solver context
718: - restart - integer restart value, this corresponds to the number of iterations of GMRES to perform before restarting
720: Options Database Key:
721: . -ksp_gmres_restart <positive integer> - integer restart value
723: Level: intermediate
725: Notes:
726: The default value is 30.
728: GMRES builds a Krylov subspace of increasing size, where each new vector is orthogonalized against the previous ones using a Gram-Schmidt process.
729: As the size of the Krylov subspace grows, the computational cost and memory requirements increase. To mitigate this issue, GMRES methods
730: usually employ restart strategies, which involve periodically deleting the Krylov subspace and beginning to generate a new one. This can help reduce
731: the computational cost and memory usage while still maintaining convergence. The maximum size of the Krylov subspace, that is the maximum number
732: of vectors orthogonalized is called the `restart` parameter.
734: A larger restart parameter generally leads to faster convergence of GMRES but the memory usage is higher than with a smaller `restart` parameter,
735: as is the average time to perform each iteration. For more ill-conditioned problems a larger restart value may be necessary.
737: `KSPBCGS` has the advantage over `KSPGMRES` in that it does not explicitly store the Krylov space and thus does not require as much memory
738: as GMRES might need.
740: .seealso: [](ch_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESGetRestart()`,
741: `KSPFGMRES`, `KSPLGMRES`, `KSPPGMRES`, `KSPAGMRES`, `KSPDGMRES`, `KSPPIPEFGMRES`
742: @*/
743: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
744: {
745: PetscFunctionBegin;
748: PetscTryMethod(ksp, "KSPGMRESSetRestart_C", (KSP, PetscInt), (ksp, restart));
749: PetscFunctionReturn(PETSC_SUCCESS);
750: }
752: /*@
753: KSPGMRESGetRestart - Gets number of iterations at which GMRES (`KSPGMRES`, `KSPFGMRES`, `KSPPGMRES`, `KSPAGMRES`, `KSPDGMRES`, `KSPPIPEFGMRES`,
754: and `KSPLGMRES`) restarts.
756: Not Collective
758: Input Parameter:
759: . ksp - the Krylov space solver context
761: Output Parameter:
762: . restart - integer restart value
764: Level: intermediate
766: .seealso: [](ch_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetRestart()`,
767: `KSPFGMRES`, `KSPLGMRES`, `KSPPGMRES`, `KSPAGMRES`, `KSPDGMRES`, `KSPPIPEFGMRES`
768: @*/
769: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
770: {
771: PetscFunctionBegin;
772: PetscUseMethod(ksp, "KSPGMRESGetRestart_C", (KSP, PetscInt *), (ksp, restart));
773: PetscFunctionReturn(PETSC_SUCCESS);
774: }
776: /*@
777: KSPGMRESSetHapTol - Sets the tolerance for detecting a happy ending in GMRES (`KSPGMRES`, `KSPFGMRES` and `KSPLGMRES` and others)
779: Logically Collective
781: Input Parameters:
782: + ksp - the Krylov space solver context
783: - tol - the tolerance for detecting a happy ending
785: Options Database Key:
786: . -ksp_gmres_haptol <positive real value> - set tolerance for determining happy breakdown
788: Level: intermediate
790: Note:
791: Happy ending is the rare case in `KSPGMRES` where a very near zero matrix entry is generated in the upper Hessenberg matrix indicating
792: an 'exact' solution has been obtained. If you attempt more iterations after this point with GMRES unstable
793: things can happen.
795: The default tolerance value for detecting a happy ending with GMRES in PETSc is 1.0e-30.
797: .seealso: [](ch_ksp), `KSPGMRES`, `KSPSetTolerances()`
798: @*/
799: PetscErrorCode KSPGMRESSetHapTol(KSP ksp, PetscReal tol)
800: {
801: PetscFunctionBegin;
803: PetscTryMethod(ksp, "KSPGMRESSetHapTol_C", (KSP, PetscReal), (ksp, tol));
804: PetscFunctionReturn(PETSC_SUCCESS);
805: }
807: /*@
808: KSPGMRESSetBreakdownTolerance - Sets the tolerance for determining divergence breakdown in `KSPGMRES` at restart.
810: Logically Collective
812: Input Parameters:
813: + ksp - the Krylov space solver context
814: - tol - the tolerance
816: Options Database Key:
817: . -ksp_gmres_breakdown_tolerance <positive real value> - set tolerance for determining divergence breakdown
819: Level: intermediate
821: Note:
822: Divergence breakdown occurs when the norm of the GMRES residual increases significantly at a restart.
823: This is defined to be $ | truenorm - gmresnorm | > tol * gmresnorm $ where $ gmresnorm $ is the norm computed
824: by the GMRES process at a restart iteration using the standard GMRES recursion formula and $ truenorm $ is computed after
825: the restart using the definition $ \| r \| = \| b - A x \|$.
827: Divergence breakdown stops the iterative solve with a `KSPConvergedReason` of `KSP_DIVERGED_BREAKDOWN` indicating the
828: GMRES solver has not converged.
830: Divergence breakdown can occur when there is an error (bug) in either the application of the matrix or the preconditioner,
831: or the preconditioner is extremely ill-conditioned.
833: The default is .1
835: .seealso: [](ch_ksp), `KSPGMRES`, `KSPSetTolerances()`, `KSPGMRESSetHapTol()`, `KSPConvergedReason`
836: @*/
837: PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP ksp, PetscReal tol)
838: {
839: PetscFunctionBegin;
841: PetscTryMethod(ksp, "KSPGMRESSetBreakdownTolerance_C", (KSP, PetscReal), (ksp, tol));
842: PetscFunctionReturn(PETSC_SUCCESS);
843: }
845: /*MC
846: KSPGMRES - Implements the Generalized Minimal Residual method {cite}`saad.schultz:gmres` with restart for solving linear systems using `KSP`.
848: Options Database Keys:
849: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
850: . -ksp_gmres_haptol <tol> - sets the tolerance for happy ending (exact convergence) of `KSPGMRES`
851: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
852: vectors are allocated as needed), see `KSPGMRESSetPreAllocateVectors()`
853: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against
854: the Krylov space (fast) (the default)
855: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
856: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
857: stability of the classical Gram-Schmidt orthogonalization.
858: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
860: Level: beginner
862: Notes:
863: Left and right preconditioning are supported, but not symmetric preconditioning.
865: Using `KSPGMRESSetPreAllocateVectors()` or `-ksp_gmres_preallocate` can improve the efficiency of the orthogonalization step with certain vector implementations.
867: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`, `KSPPGMRES`, `KSPAGMRES`, `KSPDGMRES`, `KSPPIPEFGMRES`,
868: `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
869: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
870: `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`
871: M*/
873: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
874: {
875: KSP_GMRES *gmres;
877: PetscFunctionBegin;
878: PetscCall(PetscNew(&gmres));
879: ksp->data = (void *)gmres;
881: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 4));
882: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 3));
883: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_SYMMETRIC, 2));
884: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));
885: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
887: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
888: ksp->ops->setup = KSPSetUp_GMRES;
889: ksp->ops->solve = KSPSolve_GMRES;
890: ksp->ops->reset = KSPReset_GMRES;
891: ksp->ops->destroy = KSPDestroy_GMRES;
892: ksp->ops->view = KSPView_GMRES;
893: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
894: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
895: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
896: ksp->ops->computeritz = KSPComputeRitz_GMRES;
897: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES));
898: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES));
899: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", KSPGMRESGetOrthogonalization_GMRES));
900: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES));
901: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_GMRES));
902: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES));
903: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetBreakdownTolerance_C", KSPGMRESSetBreakdownTolerance_GMRES));
904: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES));
905: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", KSPGMRESGetCGSRefinementType_GMRES));
907: gmres->haptol = 1.0e-30;
908: gmres->breakdowntol = 0.1;
909: gmres->q_preallocate = 0;
910: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
911: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
912: gmres->nrs = NULL;
913: gmres->sol_temp = NULL;
914: gmres->max_k = GMRES_DEFAULT_MAXK;
915: gmres->Rsvd = NULL;
916: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
917: gmres->orthogwork = NULL;
918: PetscFunctionReturn(PETSC_SUCCESS);
919: }