Actual source code: gmres.c
petsc-3.14.1 2020-11-03
2: /*
3: This file implements GMRES (a Generalized Minimal Residual) method.
4: Reference: Saad and Schultz, 1986.
7: Some comments on left vs. right preconditioning, and restarts.
8: Left and right preconditioning.
9: If right preconditioning is chosen, then the problem being solved
10: by gmres is actually
11: My = AB^-1 y = f
12: so the initial residual is
13: r = f - Mx
14: Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
15: residual is
16: r = f - A x
17: The final solution is then
18: x = B^-1 y
20: If left preconditioning is chosen, then the problem being solved is
21: My = B^-1 A x = B^-1 f,
22: and the initial residual is
23: r = B^-1(f - Ax)
25: Restarts: Restarts are basically solves with x0 not equal to zero.
26: Note that we can eliminate an extra application of B^-1 between
27: restarts as long as we don't require that the solution at the end
28: of an unsuccessful gmres iteration always be the solution x.
29: */
31: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
32: #define GMRES_DELTA_DIRECTIONS 10
33: #define GMRES_DEFAULT_MAXK 30
34: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP,PetscInt,PetscBool,PetscReal*);
35: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
37: PetscErrorCode KSPSetUp_GMRES(KSP ksp)
38: {
39: PetscInt hh,hes,rs,cc;
41: PetscInt max_k,k;
42: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
45: max_k = gmres->max_k; /* restart size */
46: hh = (max_k + 2) * (max_k + 1);
47: hes = (max_k + 1) * (max_k + 1);
48: rs = (max_k + 2);
49: cc = (max_k + 1);
51: PetscCalloc5(hh,&gmres->hh_origin,hes,&gmres->hes_origin,rs,&gmres->rs_origin,cc,&gmres->cc_origin,cc,&gmres->ss_origin);
52: PetscLogObjectMemory((PetscObject)ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));
54: if (ksp->calc_sings) {
55: /* Allocate workspace to hold Hessenberg matrix needed by lapack */
56: PetscMalloc1((max_k + 3)*(max_k + 9),&gmres->Rsvd);
57: PetscLogObjectMemory((PetscObject)ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
58: PetscMalloc1(6*(max_k+2),&gmres->Dsvd);
59: PetscLogObjectMemory((PetscObject)ksp,6*(max_k+2)*sizeof(PetscReal));
60: }
62: /* Allocate array to hold pointers to user vectors. Note that we need
63: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
64: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k + gmres->nextra_vecs;
66: PetscMalloc1(gmres->vecs_allocated,&gmres->vecs);
67: PetscMalloc1(VEC_OFFSET+2+max_k,&gmres->user_work);
68: PetscMalloc1(VEC_OFFSET+2+max_k,&gmres->mwork_alloc);
69: PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+2+max_k)*(sizeof(Vec*)+sizeof(PetscInt)) + gmres->vecs_allocated*sizeof(Vec));
71: if (gmres->q_preallocate) {
72: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
74: KSPCreateVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,NULL);
75: PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
77: gmres->mwork_alloc[0] = gmres->vv_allocated;
78: gmres->nwork_alloc = 1;
79: for (k=0; k<gmres->vv_allocated; k++) {
80: gmres->vecs[k] = gmres->user_work[0][k];
81: }
82: } else {
83: gmres->vv_allocated = 5;
85: KSPCreateVecs(ksp,5,&gmres->user_work[0],0,NULL);
86: PetscLogObjectParents(ksp,5,gmres->user_work[0]);
88: gmres->mwork_alloc[0] = 5;
89: gmres->nwork_alloc = 1;
90: for (k=0; k<gmres->vv_allocated; k++) {
91: gmres->vecs[k] = gmres->user_work[0][k];
92: }
93: }
94: return(0);
95: }
97: /*
98: Run gmres, possibly with restart. Return residual history if requested.
99: input parameters:
101: . gmres - structure containing parameters and work areas
103: output parameters:
104: . nres - residuals (from preconditioned system) at each step.
105: If restarting, consider passing nres+it. If null,
106: ignored
107: . itcount - number of iterations used. nres[0] to nres[itcount]
108: are defined. If null, ignored.
110: Notes:
111: On entry, the value in vector VEC_VV(0) should be the initial residual
112: (this allows shortcuts where the initial preconditioned residual is 0).
113: */
114: PetscErrorCode KSPGMRESCycle(PetscInt *itcount,KSP ksp)
115: {
116: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
117: PetscReal res,hapbnd,tt;
119: PetscInt it = 0, max_k = gmres->max_k;
120: PetscBool hapend = PETSC_FALSE;
123: if (itcount) *itcount = 0;
124: VecNormalize(VEC_VV(0),&res);
125: KSPCheckNorm(ksp,res);
127: /* the constant .1 is arbitrary, just some measure at how incorrect the residuals are */
128: if ((ksp->rnorm > 0.0) && (PetscAbsReal(res-ksp->rnorm) > .1*gmres->rnorm0)) {
129: if (ksp->errorifnotconverged) SETERRQ3(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",(double)ksp->rnorm,(double)res,(double)gmres->rnorm0);
130: else {
131: PetscInfo3(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);
132: ksp->reason = KSP_DIVERGED_BREAKDOWN;
133: return(0);
134: }
135: }
136: *GRS(0) = gmres->rnorm0 = res;
138: /* check for the convergence */
139: PetscObjectSAWsTakeAccess((PetscObject)ksp);
140: ksp->rnorm = res;
141: PetscObjectSAWsGrantAccess((PetscObject)ksp);
142: gmres->it = (it - 1);
143: KSPLogResidualHistory(ksp,res);
144: KSPMonitor(ksp,ksp->its,res);
145: if (!res) {
146: ksp->reason = KSP_CONVERGED_ATOL;
147: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
148: return(0);
149: }
151: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
152: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
153: if (it) {
154: KSPLogResidualHistory(ksp,res);
155: KSPMonitor(ksp,ksp->its,res);
156: }
157: gmres->it = (it - 1);
158: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
159: KSPGMRESGetNewVectors(ksp,it+1);
160: }
161: KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
163: /* update hessenberg matrix and do Gram-Schmidt */
164: (*gmres->orthog)(ksp,it);
165: if (ksp->reason) break;
167: /* vv(i+1) . vv(i+1) */
168: VecNormalize(VEC_VV(it+1),&tt);
169: KSPCheckNorm(ksp,tt);
171: /* save the magnitude */
172: *HH(it+1,it) = tt;
173: *HES(it+1,it) = tt;
175: /* check for the happy breakdown */
176: hapbnd = PetscAbsScalar(tt / *GRS(it));
177: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
178: if (tt < hapbnd) {
179: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %14.12e tt = %14.12e\n",(double)hapbnd,(double)tt);
180: hapend = PETSC_TRUE;
181: }
182: KSPGMRESUpdateHessenberg(ksp,it,hapend,&res);
184: it++;
185: gmres->it = (it-1); /* For converged */
186: ksp->its++;
187: ksp->rnorm = res;
188: if (ksp->reason) break;
190: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
192: /* Catch error in happy breakdown and signal convergence and break from loop */
193: if (hapend) {
194: if (ksp->normtype == KSP_NORM_NONE) { /* convergence test was skipped in this case */
195: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
196: } else if (!ksp->reason) {
197: if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",(double)res);
198: else {
199: ksp->reason = KSP_DIVERGED_BREAKDOWN;
200: break;
201: }
202: }
203: }
204: }
206: /* Monitor if we know that we will not return for a restart */
207: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
208: KSPLogResidualHistory(ksp,res);
209: KSPMonitor(ksp,ksp->its,res);
210: }
212: if (itcount) *itcount = it;
215: /*
216: Down here we have to solve for the "best" coefficients of the Krylov
217: columns, add the solution values together, and possibly unwind the
218: preconditioning from the solution
219: */
220: /* Form the solution (or the solution so far) */
221: KSPGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
222: return(0);
223: }
225: PetscErrorCode KSPSolve_GMRES(KSP ksp)
226: {
228: PetscInt its,itcount,i;
229: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
230: PetscBool guess_zero = ksp->guess_zero;
231: PetscInt N = gmres->max_k + 1;
234: if (ksp->calc_sings && !gmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
236: PetscObjectSAWsTakeAccess((PetscObject)ksp);
237: ksp->its = 0;
238: PetscObjectSAWsGrantAccess((PetscObject)ksp);
240: itcount = 0;
241: gmres->fullcycle = 0;
242: ksp->reason = KSP_CONVERGED_ITERATING;
243: ksp->rnorm = -1.0; /* special marker for KSPGMRESCycle() */
244: while (!ksp->reason) {
245: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
246: KSPGMRESCycle(&its,ksp);
247: /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
248: if the cycle is complete for the computation of the Ritz pairs */
249: if (its == gmres->max_k) {
250: gmres->fullcycle++;
251: if (ksp->calc_ritz) {
252: if (!gmres->hes_ritz) {
253: PetscMalloc1(N*N,&gmres->hes_ritz);
254: PetscLogObjectMemory((PetscObject)ksp,N*N*sizeof(PetscScalar));
255: VecDuplicateVecs(VEC_VV(0),N,&gmres->vecb);
256: }
257: PetscArraycpy(gmres->hes_ritz,gmres->hes_origin,N*N);
258: for (i=0; i<gmres->max_k+1; i++) {
259: VecCopy(VEC_VV(i),gmres->vecb[i]);
260: }
261: }
262: }
263: itcount += its;
264: if (itcount >= ksp->max_it) {
265: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
266: break;
267: }
268: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
269: }
270: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
271: return(0);
272: }
274: PetscErrorCode KSPReset_GMRES(KSP ksp)
275: {
276: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
278: PetscInt i;
281: /* Free the Hessenberg matrices */
282: PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin);
283: PetscFree(gmres->hes_ritz);
285: /* free work vectors */
286: PetscFree(gmres->vecs);
287: for (i=0; i<gmres->nwork_alloc; i++) {
288: VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
289: }
290: gmres->nwork_alloc = 0;
291: if (gmres->vecb) {
292: VecDestroyVecs(gmres->max_k+1,&gmres->vecb);
293: }
295: PetscFree(gmres->user_work);
296: PetscFree(gmres->mwork_alloc);
297: PetscFree(gmres->nrs);
298: VecDestroy(&gmres->sol_temp);
299: PetscFree(gmres->Rsvd);
300: PetscFree(gmres->Dsvd);
301: PetscFree(gmres->orthogwork);
303: gmres->vv_allocated = 0;
304: gmres->vecs_allocated = 0;
305: gmres->sol_temp = NULL;
306: return(0);
307: }
309: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
310: {
314: KSPReset_GMRES(ksp);
315: PetscFree(ksp->data);
316: /* clear composed functions */
317: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",NULL);
318: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",NULL);
319: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",NULL);
320: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",NULL);
321: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",NULL);
322: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",NULL);
323: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",NULL);
324: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",NULL);
325: return(0);
326: }
327: /*
328: KSPGMRESBuildSoln - create the solution from the starting vector and the
329: current iterates.
331: Input parameters:
332: nrs - work area of size it + 1.
333: vs - index of initial guess
334: vdest - index of result. Note that vs may == vdest (replace
335: guess with the solution).
337: This is an internal routine that knows about the GMRES internals.
338: */
339: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
340: {
341: PetscScalar tt;
343: PetscInt ii,k,j;
344: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
347: /* Solve for solution vector that minimizes the residual */
349: /* If it is < 0, no gmres steps have been performed */
350: if (it < 0) {
351: VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
352: return(0);
353: }
354: if (*HH(it,it) != 0.0) {
355: nrs[it] = *GRS(it) / *HH(it,it);
356: } else {
357: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the break down in GMRES; HH(it,it) = 0");
358: else ksp->reason = KSP_DIVERGED_BREAKDOWN;
360: PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D GRS(it) = %g\n",it,(double)PetscAbsScalar(*GRS(it)));
361: return(0);
362: }
363: for (ii=1; ii<=it; ii++) {
364: k = it - ii;
365: tt = *GRS(k);
366: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
367: if (*HH(k,k) == 0.0) {
368: if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
369: else {
370: ksp->reason = KSP_DIVERGED_BREAKDOWN;
371: PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
372: return(0);
373: }
374: }
375: nrs[k] = tt / *HH(k,k);
376: }
378: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
379: VecSet(VEC_TEMP,0.0);
380: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
382: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
383: /* add solution to previous solution */
384: if (vdest != vs) {
385: VecCopy(vs,vdest);
386: }
387: VecAXPY(vdest,1.0,VEC_TEMP);
388: return(0);
389: }
390: /*
391: Do the scalar work for the orthogonalization. Return new residual norm.
392: */
393: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
394: {
395: PetscScalar *hh,*cc,*ss,tt;
396: PetscInt j;
397: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
400: hh = HH(0,it);
401: cc = CC(0);
402: ss = SS(0);
404: /* Apply all the previously computed plane rotations to the new column
405: of the Hessenberg matrix */
406: for (j=1; j<=it; j++) {
407: tt = *hh;
408: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
409: hh++;
410: *hh = *cc++ * *hh - (*ss++ * tt);
411: }
413: /*
414: compute the new plane rotation, and apply it to:
415: 1) the right-hand-side of the Hessenberg system
416: 2) the new column of the Hessenberg matrix
417: thus obtaining the updated value of the residual
418: */
419: if (!hapend) {
420: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
421: if (tt == 0.0) {
422: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"tt == 0.0");
423: else {
424: ksp->reason = KSP_DIVERGED_NULL;
425: return(0);
426: }
427: }
428: *cc = *hh / tt;
429: *ss = *(hh+1) / tt;
430: *GRS(it+1) = -(*ss * *GRS(it));
431: *GRS(it) = PetscConj(*cc) * *GRS(it);
432: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
433: *res = PetscAbsScalar(*GRS(it+1));
434: } else {
435: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
436: another rotation matrix (so RH doesn't change). The new residual is
437: always the new sine term times the residual from last time (GRS(it)),
438: but now the new sine rotation would be zero...so the residual should
439: be zero...so we will multiply "zero" by the last residual. This might
440: not be exactly what we want to do here -could just return "zero". */
442: *res = 0.0;
443: }
444: return(0);
445: }
446: /*
447: This routine allocates more work vectors, starting from VEC_VV(it).
448: */
449: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
450: {
451: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
453: PetscInt nwork = gmres->nwork_alloc,k,nalloc;
456: nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
457: /* Adjust the number to allocate to make sure that we don't exceed the
458: number of available slots */
459: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) {
460: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
461: }
462: if (!nalloc) return(0);
464: gmres->vv_allocated += nalloc;
466: KSPCreateVecs(ksp,nalloc,&gmres->user_work[nwork],0,NULL);
467: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
469: gmres->mwork_alloc[nwork] = nalloc;
470: for (k=0; k<nalloc; k++) {
471: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
472: }
473: gmres->nwork_alloc++;
474: return(0);
475: }
477: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
478: {
479: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
483: if (!ptr) {
484: if (!gmres->sol_temp) {
485: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
486: PetscLogObjectParent((PetscObject)ksp,(PetscObject)gmres->sol_temp);
487: }
488: ptr = gmres->sol_temp;
489: }
490: if (!gmres->nrs) {
491: /* allocate the work area */
492: PetscMalloc1(gmres->max_k,&gmres->nrs);
493: PetscLogObjectMemory((PetscObject)ksp,gmres->max_k);
494: }
496: KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
497: if (result) *result = ptr;
498: return(0);
499: }
501: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
502: {
503: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
504: const char *cstr;
506: PetscBool iascii,isstring;
509: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
510: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
511: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
512: switch (gmres->cgstype) {
513: case (KSP_GMRES_CGS_REFINE_NEVER):
514: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
515: break;
516: case (KSP_GMRES_CGS_REFINE_ALWAYS):
517: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
518: break;
519: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
520: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
521: break;
522: default:
523: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
524: }
525: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
526: cstr = "Modified Gram-Schmidt Orthogonalization";
527: } else {
528: cstr = "unknown orthogonalization";
529: }
530: if (iascii) {
531: PetscViewerASCIIPrintf(viewer," restart=%D, using %s\n",gmres->max_k,cstr);
532: PetscViewerASCIIPrintf(viewer," happy breakdown tolerance %g\n",(double)gmres->haptol);
533: } else if (isstring) {
534: PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
535: }
536: return(0);
537: }
539: /*@C
540: KSPGMRESMonitorKrylov - Calls VecView() for each new direction in the GMRES accumulated Krylov space.
542: Collective on ksp
544: Input Parameters:
545: + ksp - the KSP context
546: . its - iteration number
547: . fgnorm - 2-norm of residual (or gradient)
548: - dummy - an collection of viewers created with KSPViewerCreate()
550: Options Database Keys:
551: . -ksp_gmres_kyrlov_monitor
553: Notes:
554: A new PETSCVIEWERDRAW is created for each Krylov vector so they can all be simultaneously viewed
555: Level: intermediate
557: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), KSPViewersCreate(), KSPViewersDestroy()
558: @*/
559: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
560: {
561: PetscViewers viewers = (PetscViewers)dummy;
562: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
564: Vec x;
565: PetscViewer viewer;
566: PetscBool flg;
569: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
570: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
571: if (!flg) {
572: PetscViewerSetType(viewer,PETSCVIEWERDRAW);
573: PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
574: }
575: x = VEC_VV(gmres->it+1);
576: VecView(x,viewer);
577: return(0);
578: }
580: PetscErrorCode KSPSetFromOptions_GMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
581: {
583: PetscInt restart;
584: PetscReal haptol;
585: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
586: PetscBool flg;
589: PetscOptionsHead(PetscOptionsObject,"KSP GMRES Options");
590: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
591: if (flg) { KSPGMRESSetRestart(ksp,restart); }
592: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
593: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
594: flg = PETSC_FALSE;
595: PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,NULL);
596: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
597: PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
598: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
599: PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
600: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
601: PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
602: KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
603: flg = PETSC_FALSE;
604: PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,NULL);
605: if (flg) {
606: PetscViewers viewers;
607: PetscViewersCreate(PetscObjectComm((PetscObject)ksp),&viewers);
608: KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
609: }
610: PetscOptionsTail();
611: return(0);
612: }
614: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
615: {
616: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
619: if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
620: gmres->haptol = tol;
621: return(0);
622: }
624: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
625: {
626: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
629: *max_k = gmres->max_k;
630: return(0);
631: }
633: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
634: {
635: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
639: if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
640: if (!ksp->setupstage) {
641: gmres->max_k = max_k;
642: } else if (gmres->max_k != max_k) {
643: gmres->max_k = max_k;
644: ksp->setupstage = KSP_SETUP_NEW;
645: /* free the data structures, then create them again */
646: KSPReset_GMRES(ksp);
647: }
648: return(0);
649: }
651: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
652: {
654: ((KSP_GMRES*)ksp->data)->orthog = fcn;
655: return(0);
656: }
658: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
659: {
661: *fcn = ((KSP_GMRES*)ksp->data)->orthog;
662: return(0);
663: }
665: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
666: {
667: KSP_GMRES *gmres;
670: gmres = (KSP_GMRES*)ksp->data;
671: gmres->q_preallocate = 1;
672: return(0);
673: }
675: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
676: {
677: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
680: gmres->cgstype = type;
681: return(0);
682: }
684: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
685: {
686: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
689: *type = gmres->cgstype;
690: return(0);
691: }
693: /*@
694: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
695: in the classical Gram Schmidt orthogonalization.
697: Logically Collective on ksp
699: Input Parameters:
700: + ksp - the Krylov space context
701: - type - the type of refinement
703: Options Database:
704: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
706: Level: intermediate
708: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
709: KSPGMRESGetOrthogonalization()
710: @*/
711: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
712: {
718: PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
719: return(0);
720: }
722: /*@
723: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
724: in the classical Gram Schmidt orthogonalization.
726: Not Collective
728: Input Parameter:
729: . ksp - the Krylov space context
731: Output Parameter:
732: . type - the type of refinement
734: Options Database:
735: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
737: Level: intermediate
739: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
740: KSPGMRESGetOrthogonalization()
741: @*/
742: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
743: {
748: PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
749: return(0);
750: }
753: /*@
754: KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.
756: Logically Collective on ksp
758: Input Parameters:
759: + ksp - the Krylov space context
760: - restart - integer restart value
762: Options Database:
763: . -ksp_gmres_restart <positive integer>
765: Note: The default value is 30.
767: Level: intermediate
769: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
770: @*/
771: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
772: {
778: PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
779: return(0);
780: }
782: /*@
783: KSPGMRESGetRestart - Gets number of iterations at which GMRES, FGMRES and LGMRES restarts.
785: Not Collective
787: Input Parameter:
788: . ksp - the Krylov space context
790: Output Parameter:
791: . restart - integer restart value
793: Note: The default value is 30.
795: Level: intermediate
797: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
798: @*/
799: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
800: {
804: PetscUseMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
805: return(0);
806: }
808: /*@
809: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
811: Logically Collective on ksp
813: Input Parameters:
814: + ksp - the Krylov space context
815: - tol - the tolerance
817: Options Database:
818: . -ksp_gmres_haptol <positive real value>
820: Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
821: a certain number of iterations. If you attempt more iterations after this point unstable
822: things can happen hence very occasionally you may need to set this value to detect this condition
824: Level: intermediate
826: .seealso: KSPSetTolerances()
827: @*/
828: PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
829: {
834: PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
835: return(0);
836: }
838: /*MC
839: KSPGMRES - Implements the Generalized Minimal Residual method.
840: (Saad and Schultz, 1986) with restart
843: Options Database Keys:
844: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
845: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
846: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
847: vectors are allocated as needed)
848: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
849: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
850: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
851: stability of the classical Gram-Schmidt orthogonalization.
852: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
854: Level: beginner
856: Notes:
857: Left and right preconditioning are supported, but not symmetric preconditioning.
859: References:
860: . 1. - YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
861: SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.
863: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
864: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
865: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
866: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()
868: M*/
870: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
871: {
872: KSP_GMRES *gmres;
876: PetscNewLog(ksp,&gmres);
877: ksp->data = (void*)gmres;
879: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,4);
880: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
881: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,2);
882: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
883: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
885: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
886: ksp->ops->setup = KSPSetUp_GMRES;
887: ksp->ops->solve = KSPSolve_GMRES;
888: ksp->ops->reset = KSPReset_GMRES;
889: ksp->ops->destroy = KSPDestroy_GMRES;
890: ksp->ops->view = KSPView_GMRES;
891: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
892: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
893: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
894: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_ESSL)
895: ksp->ops->computeritz = KSPComputeRitz_GMRES;
896: #endif
897: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
898: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
899: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
900: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
901: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
902: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
903: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
904: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
906: gmres->haptol = 1.0e-30;
907: gmres->q_preallocate = 0;
908: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
909: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
910: gmres->nrs = NULL;
911: gmres->sol_temp = NULL;
912: gmres->max_k = GMRES_DEFAULT_MAXK;
913: gmres->Rsvd = NULL;
914: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
915: gmres->orthogwork = NULL;
916: return(0);
917: }