Actual source code: gmres.c
petsc-3.11.4 2019-09-28
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_norm,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_norm);
125: KSPCheckNorm(ksp,res_norm);
126: res = res_norm;
127: *GRS(0) = res_norm;
129: /* check for the convergence */
130: PetscObjectSAWsTakeAccess((PetscObject)ksp);
131: ksp->rnorm = res;
132: PetscObjectSAWsGrantAccess((PetscObject)ksp);
133: gmres->it = (it - 1);
134: KSPLogResidualHistory(ksp,res);
135: KSPMonitor(ksp,ksp->its,res);
136: if (!res) {
137: ksp->reason = KSP_CONVERGED_ATOL;
138: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
139: return(0);
140: }
142: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
143: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
144: if (it) {
145: KSPLogResidualHistory(ksp,res);
146: KSPMonitor(ksp,ksp->its,res);
147: }
148: gmres->it = (it - 1);
149: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
150: KSPGMRESGetNewVectors(ksp,it+1);
151: }
152: KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
154: /* update hessenberg matrix and do Gram-Schmidt */
155: (*gmres->orthog)(ksp,it);
156: if (ksp->reason) break;
158: /* vv(i+1) . vv(i+1) */
159: VecNormalize(VEC_VV(it+1),&tt);
160: KSPCheckNorm(ksp,tt);
162: /* save the magnitude */
163: *HH(it+1,it) = tt;
164: *HES(it+1,it) = tt;
166: /* check for the happy breakdown */
167: hapbnd = PetscAbsScalar(tt / *GRS(it));
168: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
169: if (tt < hapbnd) {
170: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %14.12e tt = %14.12e\n",(double)hapbnd,(double)tt);
171: hapend = PETSC_TRUE;
172: }
173: KSPGMRESUpdateHessenberg(ksp,it,hapend,&res);
175: it++;
176: gmres->it = (it-1); /* For converged */
177: ksp->its++;
178: ksp->rnorm = res;
179: if (ksp->reason) break;
181: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
183: /* Catch error in happy breakdown and signal convergence and break from loop */
184: if (hapend) {
185: if (ksp->normtype == KSP_NORM_NONE) { /* convergence test was skipped in this case */
186: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
187: } else if (!ksp->reason) {
188: 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);
189: else {
190: ksp->reason = KSP_DIVERGED_BREAKDOWN;
191: break;
192: }
193: }
194: }
195: }
197: /* Monitor if we know that we will not return for a restart */
198: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
199: KSPLogResidualHistory(ksp,res);
200: KSPMonitor(ksp,ksp->its,res);
201: }
203: if (itcount) *itcount = it;
206: /*
207: Down here we have to solve for the "best" coefficients of the Krylov
208: columns, add the solution values together, and possibly unwind the
209: preconditioning from the solution
210: */
211: /* Form the solution (or the solution so far) */
212: KSPGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
213: return(0);
214: }
216: PetscErrorCode KSPSolve_GMRES(KSP ksp)
217: {
219: PetscInt its,itcount,i;
220: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
221: PetscBool guess_zero = ksp->guess_zero;
222: PetscInt N = gmres->max_k + 1;
223: PetscBLASInt bN;
226: if (ksp->calc_sings && !gmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
228: PetscObjectSAWsTakeAccess((PetscObject)ksp);
229: ksp->its = 0;
230: PetscObjectSAWsGrantAccess((PetscObject)ksp);
232: itcount = 0;
233: gmres->fullcycle = 0;
234: ksp->reason = KSP_CONVERGED_ITERATING;
235: while (!ksp->reason) {
236: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
237: KSPGMRESCycle(&its,ksp);
238: /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
239: if the cycle is complete for the computation of the Ritz pairs */
240: if (its == gmres->max_k) {
241: gmres->fullcycle++;
242: if (ksp->calc_ritz) {
243: if (!gmres->hes_ritz) {
244: PetscMalloc1(N*N,&gmres->hes_ritz);
245: PetscLogObjectMemory((PetscObject)ksp,N*N*sizeof(PetscScalar));
246: VecDuplicateVecs(VEC_VV(0),N,&gmres->vecb);
247: }
248: PetscBLASIntCast(N,&bN);
249: PetscMemcpy(gmres->hes_ritz,gmres->hes_origin,bN*bN*sizeof(PetscReal));
250: for (i=0; i<gmres->max_k+1; i++) {
251: VecCopy(VEC_VV(i),gmres->vecb[i]);
252: }
253: }
254: }
255: itcount += its;
256: if (itcount >= ksp->max_it) {
257: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
258: break;
259: }
260: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
261: }
262: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
263: return(0);
264: }
266: PetscErrorCode KSPReset_GMRES(KSP ksp)
267: {
268: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
270: PetscInt i;
273: /* Free the Hessenberg matrices */
274: PetscFree6(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin,gmres->hes_ritz);
276: /* free work vectors */
277: PetscFree(gmres->vecs);
278: for (i=0; i<gmres->nwork_alloc; i++) {
279: VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
280: }
281: gmres->nwork_alloc = 0;
282: if (gmres->vecb) {
283: VecDestroyVecs(gmres->max_k+1,&gmres->vecb);
284: }
286: PetscFree(gmres->user_work);
287: PetscFree(gmres->mwork_alloc);
288: PetscFree(gmres->nrs);
289: VecDestroy(&gmres->sol_temp);
290: PetscFree(gmres->Rsvd);
291: PetscFree(gmres->Dsvd);
292: PetscFree(gmres->orthogwork);
294: gmres->sol_temp = 0;
295: gmres->vv_allocated = 0;
296: gmres->vecs_allocated = 0;
297: gmres->sol_temp = 0;
298: return(0);
299: }
301: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
302: {
306: KSPReset_GMRES(ksp);
307: PetscFree(ksp->data);
308: /* clear composed functions */
309: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",NULL);
310: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",NULL);
311: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",NULL);
312: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",NULL);
313: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",NULL);
314: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",NULL);
315: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",NULL);
316: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",NULL);
317: return(0);
318: }
319: /*
320: KSPGMRESBuildSoln - create the solution from the starting vector and the
321: current iterates.
323: Input parameters:
324: nrs - work area of size it + 1.
325: vs - index of initial guess
326: vdest - index of result. Note that vs may == vdest (replace
327: guess with the solution).
329: This is an internal routine that knows about the GMRES internals.
330: */
331: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
332: {
333: PetscScalar tt;
335: PetscInt ii,k,j;
336: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
339: /* Solve for solution vector that minimizes the residual */
341: /* If it is < 0, no gmres steps have been performed */
342: if (it < 0) {
343: VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
344: return(0);
345: }
346: if (*HH(it,it) != 0.0) {
347: nrs[it] = *GRS(it) / *HH(it,it);
348: } else {
349: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the break down in GMRES; HH(it,it) = 0");
350: else ksp->reason = KSP_DIVERGED_BREAKDOWN;
352: 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)));
353: return(0);
354: }
355: for (ii=1; ii<=it; ii++) {
356: k = it - ii;
357: tt = *GRS(k);
358: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
359: if (*HH(k,k) == 0.0) {
360: 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);
361: else {
362: ksp->reason = KSP_DIVERGED_BREAKDOWN;
363: PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
364: return(0);
365: }
366: }
367: nrs[k] = tt / *HH(k,k);
368: }
370: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
371: VecSet(VEC_TEMP,0.0);
372: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
374: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
375: /* add solution to previous solution */
376: if (vdest != vs) {
377: VecCopy(vs,vdest);
378: }
379: VecAXPY(vdest,1.0,VEC_TEMP);
380: return(0);
381: }
382: /*
383: Do the scalar work for the orthogonalization. Return new residual norm.
384: */
385: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
386: {
387: PetscScalar *hh,*cc,*ss,tt;
388: PetscInt j;
389: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
392: hh = HH(0,it);
393: cc = CC(0);
394: ss = SS(0);
396: /* Apply all the previously computed plane rotations to the new column
397: of the Hessenberg matrix */
398: for (j=1; j<=it; j++) {
399: tt = *hh;
400: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
401: hh++;
402: *hh = *cc++ * *hh - (*ss++ * tt);
403: }
405: /*
406: compute the new plane rotation, and apply it to:
407: 1) the right-hand-side of the Hessenberg system
408: 2) the new column of the Hessenberg matrix
409: thus obtaining the updated value of the residual
410: */
411: if (!hapend) {
412: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
413: if (tt == 0.0) {
414: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"tt == 0.0");
415: else {
416: ksp->reason = KSP_DIVERGED_NULL;
417: return(0);
418: }
419: }
420: *cc = *hh / tt;
421: *ss = *(hh+1) / tt;
422: *GRS(it+1) = -(*ss * *GRS(it));
423: *GRS(it) = PetscConj(*cc) * *GRS(it);
424: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
425: *res = PetscAbsScalar(*GRS(it+1));
426: } else {
427: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
428: another rotation matrix (so RH doesn't change). The new residual is
429: always the new sine term times the residual from last time (GRS(it)),
430: but now the new sine rotation would be zero...so the residual should
431: be zero...so we will multiply "zero" by the last residual. This might
432: not be exactly what we want to do here -could just return "zero". */
434: *res = 0.0;
435: }
436: return(0);
437: }
438: /*
439: This routine allocates more work vectors, starting from VEC_VV(it).
440: */
441: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
442: {
443: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
445: PetscInt nwork = gmres->nwork_alloc,k,nalloc;
448: nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
449: /* Adjust the number to allocate to make sure that we don't exceed the
450: number of available slots */
451: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) {
452: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
453: }
454: if (!nalloc) return(0);
456: gmres->vv_allocated += nalloc;
458: KSPCreateVecs(ksp,nalloc,&gmres->user_work[nwork],0,NULL);
459: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
461: gmres->mwork_alloc[nwork] = nalloc;
462: for (k=0; k<nalloc; k++) {
463: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
464: }
465: gmres->nwork_alloc++;
466: return(0);
467: }
469: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
470: {
471: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
475: if (!ptr) {
476: if (!gmres->sol_temp) {
477: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
478: PetscLogObjectParent((PetscObject)ksp,(PetscObject)gmres->sol_temp);
479: }
480: ptr = gmres->sol_temp;
481: }
482: if (!gmres->nrs) {
483: /* allocate the work area */
484: PetscMalloc1(gmres->max_k,&gmres->nrs);
485: PetscLogObjectMemory((PetscObject)ksp,gmres->max_k*sizeof(PetscScalar));
486: }
488: KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
489: if (result) *result = ptr;
490: return(0);
491: }
493: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
494: {
495: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
496: const char *cstr;
498: PetscBool iascii,isstring;
501: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
502: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
503: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
504: switch (gmres->cgstype) {
505: case (KSP_GMRES_CGS_REFINE_NEVER):
506: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
507: break;
508: case (KSP_GMRES_CGS_REFINE_ALWAYS):
509: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
510: break;
511: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
512: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
513: break;
514: default:
515: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
516: }
517: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
518: cstr = "Modified Gram-Schmidt Orthogonalization";
519: } else {
520: cstr = "unknown orthogonalization";
521: }
522: if (iascii) {
523: PetscViewerASCIIPrintf(viewer," restart=%D, using %s\n",gmres->max_k,cstr);
524: PetscViewerASCIIPrintf(viewer," happy breakdown tolerance %g\n",(double)gmres->haptol);
525: } else if (isstring) {
526: PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
527: }
528: return(0);
529: }
531: /*@C
532: KSPGMRESMonitorKrylov - Calls VecView() for each new direction in the GMRES accumulated Krylov space.
534: Collective on KSP
536: Input Parameters:
537: + ksp - the KSP context
538: . its - iteration number
539: . fgnorm - 2-norm of residual (or gradient)
540: - dummy - an collection of viewers created with KSPViewerCreate()
542: Options Database Keys:
543: . -ksp_gmres_kyrlov_monitor
545: Notes:
546: A new PETSCVIEWERDRAW is created for each Krylov vector so they can all be simultaneously viewed
547: Level: intermediate
549: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
551: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), KSPViewersCreate(), KSPViewersDestroy()
552: @*/
553: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
554: {
555: PetscViewers viewers = (PetscViewers)dummy;
556: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
558: Vec x;
559: PetscViewer viewer;
560: PetscBool flg;
563: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
564: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
565: if (!flg) {
566: PetscViewerSetType(viewer,PETSCVIEWERDRAW);
567: PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
568: }
569: x = VEC_VV(gmres->it+1);
570: VecView(x,viewer);
571: return(0);
572: }
574: PetscErrorCode KSPSetFromOptions_GMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
575: {
577: PetscInt restart;
578: PetscReal haptol;
579: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
580: PetscBool flg;
583: PetscOptionsHead(PetscOptionsObject,"KSP GMRES Options");
584: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
585: if (flg) { KSPGMRESSetRestart(ksp,restart); }
586: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
587: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
588: flg = PETSC_FALSE;
589: PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,NULL);
590: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
591: PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
592: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
593: PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
594: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
595: PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
596: KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
597: flg = PETSC_FALSE;
598: PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,NULL);
599: if (flg) {
600: PetscViewers viewers;
601: PetscViewersCreate(PetscObjectComm((PetscObject)ksp),&viewers);
602: KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
603: }
604: PetscOptionsTail();
605: return(0);
606: }
608: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
609: {
610: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
613: if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
614: gmres->haptol = tol;
615: return(0);
616: }
618: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
619: {
620: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
623: *max_k = gmres->max_k;
624: return(0);
625: }
627: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
628: {
629: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
633: if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
634: if (!ksp->setupstage) {
635: gmres->max_k = max_k;
636: } else if (gmres->max_k != max_k) {
637: gmres->max_k = max_k;
638: ksp->setupstage = KSP_SETUP_NEW;
639: /* free the data structures, then create them again */
640: KSPReset_GMRES(ksp);
641: }
642: return(0);
643: }
645: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
646: {
648: ((KSP_GMRES*)ksp->data)->orthog = fcn;
649: return(0);
650: }
652: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
653: {
655: *fcn = ((KSP_GMRES*)ksp->data)->orthog;
656: return(0);
657: }
659: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
660: {
661: KSP_GMRES *gmres;
664: gmres = (KSP_GMRES*)ksp->data;
665: gmres->q_preallocate = 1;
666: return(0);
667: }
669: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
670: {
671: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
674: gmres->cgstype = type;
675: return(0);
676: }
678: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
679: {
680: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
683: *type = gmres->cgstype;
684: return(0);
685: }
687: /*@
688: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
689: in the classical Gram Schmidt orthogonalization.
691: Logically Collective on KSP
693: Input Parameters:
694: + ksp - the Krylov space context
695: - type - the type of refinement
697: Options Database:
698: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
700: Level: intermediate
702: .keywords: KSP, GMRES, iterative refinement
704: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
705: KSPGMRESGetOrthogonalization()
706: @*/
707: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
708: {
714: PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
715: return(0);
716: }
718: /*@
719: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
720: in the classical Gram Schmidt orthogonalization.
722: Not Collective
724: Input Parameter:
725: . ksp - the Krylov space context
727: Output Parameter:
728: . type - the type of refinement
730: Options Database:
731: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
733: Level: intermediate
735: .keywords: KSP, GMRES, iterative refinement
737: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
738: KSPGMRESGetOrthogonalization()
739: @*/
740: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
741: {
746: PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
747: return(0);
748: }
751: /*@
752: KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.
754: Logically Collective on KSP
756: Input Parameters:
757: + ksp - the Krylov space context
758: - restart - integer restart value
760: Options Database:
761: . -ksp_gmres_restart <positive integer>
763: Note: The default value is 30.
765: Level: intermediate
767: .keywords: KSP, GMRES, restart, iterations
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: .keywords: KSP, GMRES, restart, iterations
799: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
800: @*/
801: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
802: {
806: PetscUseMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
807: return(0);
808: }
810: /*@
811: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
813: Logically Collective on KSP
815: Input Parameters:
816: + ksp - the Krylov space context
817: - tol - the tolerance
819: Options Database:
820: . -ksp_gmres_haptol <positive real value>
822: Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
823: a certain number of iterations. If you attempt more iterations after this point unstable
824: things can happen hence very occasionally you may need to set this value to detect this condition
826: Level: intermediate
828: .keywords: KSP, GMRES, tolerance
830: .seealso: KSPSetTolerances()
831: @*/
832: PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
833: {
838: PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
839: return(0);
840: }
842: /*MC
843: KSPGMRES - Implements the Generalized Minimal Residual method.
844: (Saad and Schultz, 1986) with restart
847: Options Database Keys:
848: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
849: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
850: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
851: vectors are allocated as needed)
852: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
853: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
854: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
855: stability of the classical Gram-Schmidt orthogonalization.
856: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
858: Level: beginner
860: Notes:
861: Left and right preconditioning are supported, but not symmetric preconditioning.
863: References:
864: . 1. - YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
865: SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.
867: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
868: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
869: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
870: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()
872: M*/
874: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
875: {
876: KSP_GMRES *gmres;
880: PetscNewLog(ksp,&gmres);
881: ksp->data = (void*)gmres;
883: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,4);
884: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
885: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,2);
886: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
887: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
889: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
890: ksp->ops->setup = KSPSetUp_GMRES;
891: ksp->ops->solve = KSPSolve_GMRES;
892: ksp->ops->reset = KSPReset_GMRES;
893: ksp->ops->destroy = KSPDestroy_GMRES;
894: ksp->ops->view = KSPView_GMRES;
895: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
896: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
897: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
898: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_ESSL)
899: ksp->ops->computeritz = KSPComputeRitz_GMRES;
900: #endif
901: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
902: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
903: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
904: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
905: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
906: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
907: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
908: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
910: gmres->haptol = 1.0e-30;
911: gmres->q_preallocate = 0;
912: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
913: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
914: gmres->nrs = 0;
915: gmres->sol_temp = 0;
916: gmres->max_k = GMRES_DEFAULT_MAXK;
917: gmres->Rsvd = 0;
918: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
919: gmres->orthogwork = 0;
920: return(0);
921: }