Actual source code: gmres.c
petsc-3.10.5 2019-03-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);
161: /* save the magnitude */
162: *HH(it+1,it) = tt;
163: *HES(it+1,it) = tt;
165: /* check for the happy breakdown */
166: hapbnd = PetscAbsScalar(tt / *GRS(it));
167: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
168: if (tt < hapbnd) {
169: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %14.12e tt = %14.12e\n",(double)hapbnd,(double)tt);
170: hapend = PETSC_TRUE;
171: }
172: KSPGMRESUpdateHessenberg(ksp,it,hapend,&res);
174: it++;
175: gmres->it = (it-1); /* For converged */
176: ksp->its++;
177: ksp->rnorm = res;
178: if (ksp->reason) break;
180: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
182: /* Catch error in happy breakdown and signal convergence and break from loop */
183: if (hapend) {
184: if (!ksp->reason) {
185: 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);
186: else {
187: ksp->reason = KSP_DIVERGED_BREAKDOWN;
188: break;
189: }
190: }
191: }
192: }
194: /* Monitor if we know that we will not return for a restart */
195: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
196: KSPLogResidualHistory(ksp,res);
197: KSPMonitor(ksp,ksp->its,res);
198: }
200: if (itcount) *itcount = it;
203: /*
204: Down here we have to solve for the "best" coefficients of the Krylov
205: columns, add the solution values together, and possibly unwind the
206: preconditioning from the solution
207: */
208: /* Form the solution (or the solution so far) */
209: KSPGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
210: return(0);
211: }
213: PetscErrorCode KSPSolve_GMRES(KSP ksp)
214: {
216: PetscInt its,itcount,i;
217: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
218: PetscBool guess_zero = ksp->guess_zero;
219: PetscInt N = gmres->max_k + 1;
220: PetscBLASInt bN;
223: if (ksp->calc_sings && !gmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
225: PetscObjectSAWsTakeAccess((PetscObject)ksp);
226: ksp->its = 0;
227: PetscObjectSAWsGrantAccess((PetscObject)ksp);
229: itcount = 0;
230: gmres->fullcycle = 0;
231: ksp->reason = KSP_CONVERGED_ITERATING;
232: while (!ksp->reason) {
233: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
234: KSPGMRESCycle(&its,ksp);
235: /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
236: if the cycle is complete for the computation of the Ritz pairs */
237: if (its == gmres->max_k) {
238: gmres->fullcycle++;
239: if (ksp->calc_ritz) {
240: if (!gmres->hes_ritz) {
241: PetscMalloc1(N*N,&gmres->hes_ritz);
242: PetscLogObjectMemory((PetscObject)ksp,N*N*sizeof(PetscScalar));
243: VecDuplicateVecs(VEC_VV(0),N,&gmres->vecb);
244: }
245: PetscBLASIntCast(N,&bN);
246: PetscMemcpy(gmres->hes_ritz,gmres->hes_origin,bN*bN*sizeof(PetscReal));
247: for (i=0; i<gmres->max_k+1; i++) {
248: VecCopy(VEC_VV(i),gmres->vecb[i]);
249: }
250: }
251: }
252: itcount += its;
253: if (itcount >= ksp->max_it) {
254: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
255: break;
256: }
257: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
258: }
259: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
260: return(0);
261: }
263: PetscErrorCode KSPReset_GMRES(KSP ksp)
264: {
265: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
267: PetscInt i;
270: /* Free the Hessenberg matrices */
271: PetscFree6(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin,gmres->hes_ritz);
273: /* free work vectors */
274: PetscFree(gmres->vecs);
275: for (i=0; i<gmres->nwork_alloc; i++) {
276: VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
277: }
278: gmres->nwork_alloc = 0;
279: if (gmres->vecb) {
280: VecDestroyVecs(gmres->max_k+1,&gmres->vecb);
281: }
283: PetscFree(gmres->user_work);
284: PetscFree(gmres->mwork_alloc);
285: PetscFree(gmres->nrs);
286: VecDestroy(&gmres->sol_temp);
287: PetscFree(gmres->Rsvd);
288: PetscFree(gmres->Dsvd);
289: PetscFree(gmres->orthogwork);
291: gmres->sol_temp = 0;
292: gmres->vv_allocated = 0;
293: gmres->vecs_allocated = 0;
294: gmres->sol_temp = 0;
295: return(0);
296: }
298: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
299: {
303: KSPReset_GMRES(ksp);
304: PetscFree(ksp->data);
305: /* clear composed functions */
306: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",NULL);
307: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",NULL);
308: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",NULL);
309: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",NULL);
310: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",NULL);
311: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",NULL);
312: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",NULL);
313: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",NULL);
314: return(0);
315: }
316: /*
317: KSPGMRESBuildSoln - create the solution from the starting vector and the
318: current iterates.
320: Input parameters:
321: nrs - work area of size it + 1.
322: vs - index of initial guess
323: vdest - index of result. Note that vs may == vdest (replace
324: guess with the solution).
326: This is an internal routine that knows about the GMRES internals.
327: */
328: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
329: {
330: PetscScalar tt;
332: PetscInt ii,k,j;
333: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
336: /* Solve for solution vector that minimizes the residual */
338: /* If it is < 0, no gmres steps have been performed */
339: if (it < 0) {
340: VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
341: return(0);
342: }
343: if (*HH(it,it) != 0.0) {
344: nrs[it] = *GRS(it) / *HH(it,it);
345: } else {
346: ksp->reason = KSP_DIVERGED_BREAKDOWN;
348: 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)));
349: return(0);
350: }
351: for (ii=1; ii<=it; ii++) {
352: k = it - ii;
353: tt = *GRS(k);
354: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
355: if (*HH(k,k) == 0.0) {
356: ksp->reason = KSP_DIVERGED_BREAKDOWN;
358: PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
359: return(0);
360: }
361: nrs[k] = tt / *HH(k,k);
362: }
364: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
365: VecSet(VEC_TEMP,0.0);
366: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
368: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
369: /* add solution to previous solution */
370: if (vdest != vs) {
371: VecCopy(vs,vdest);
372: }
373: VecAXPY(vdest,1.0,VEC_TEMP);
374: return(0);
375: }
376: /*
377: Do the scalar work for the orthogonalization. Return new residual norm.
378: */
379: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
380: {
381: PetscScalar *hh,*cc,*ss,tt;
382: PetscInt j;
383: KSP_GMRES *gmres = (KSP_GMRES*)(ksp->data);
386: hh = HH(0,it);
387: cc = CC(0);
388: ss = SS(0);
390: /* Apply all the previously computed plane rotations to the new column
391: of the Hessenberg matrix */
392: for (j=1; j<=it; j++) {
393: tt = *hh;
394: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
395: hh++;
396: *hh = *cc++ * *hh - (*ss++ * tt);
397: }
399: /*
400: compute the new plane rotation, and apply it to:
401: 1) the right-hand-side of the Hessenberg system
402: 2) the new column of the Hessenberg matrix
403: thus obtaining the updated value of the residual
404: */
405: if (!hapend) {
406: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
407: if (tt == 0.0) {
408: ksp->reason = KSP_DIVERGED_NULL;
409: return(0);
410: }
411: *cc = *hh / tt;
412: *ss = *(hh+1) / tt;
413: *GRS(it+1) = -(*ss * *GRS(it));
414: *GRS(it) = PetscConj(*cc) * *GRS(it);
415: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
416: *res = PetscAbsScalar(*GRS(it+1));
417: } else {
418: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
419: another rotation matrix (so RH doesn't change). The new residual is
420: always the new sine term times the residual from last time (GRS(it)),
421: but now the new sine rotation would be zero...so the residual should
422: be zero...so we will multiply "zero" by the last residual. This might
423: not be exactly what we want to do here -could just return "zero". */
425: *res = 0.0;
426: }
427: return(0);
428: }
429: /*
430: This routine allocates more work vectors, starting from VEC_VV(it).
431: */
432: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
433: {
434: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
436: PetscInt nwork = gmres->nwork_alloc,k,nalloc;
439: nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
440: /* Adjust the number to allocate to make sure that we don't exceed the
441: number of available slots */
442: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) {
443: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
444: }
445: if (!nalloc) return(0);
447: gmres->vv_allocated += nalloc;
449: KSPCreateVecs(ksp,nalloc,&gmres->user_work[nwork],0,NULL);
450: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
452: gmres->mwork_alloc[nwork] = nalloc;
453: for (k=0; k<nalloc; k++) {
454: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
455: }
456: gmres->nwork_alloc++;
457: return(0);
458: }
460: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
461: {
462: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
466: if (!ptr) {
467: if (!gmres->sol_temp) {
468: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
469: PetscLogObjectParent((PetscObject)ksp,(PetscObject)gmres->sol_temp);
470: }
471: ptr = gmres->sol_temp;
472: }
473: if (!gmres->nrs) {
474: /* allocate the work area */
475: PetscMalloc1(gmres->max_k,&gmres->nrs);
476: PetscLogObjectMemory((PetscObject)ksp,gmres->max_k*sizeof(PetscScalar));
477: }
479: KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
480: if (result) *result = ptr;
481: return(0);
482: }
484: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
485: {
486: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
487: const char *cstr;
489: PetscBool iascii,isstring;
492: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
493: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
494: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
495: switch (gmres->cgstype) {
496: case (KSP_GMRES_CGS_REFINE_NEVER):
497: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
498: break;
499: case (KSP_GMRES_CGS_REFINE_ALWAYS):
500: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
501: break;
502: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
503: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
504: break;
505: default:
506: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
507: }
508: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
509: cstr = "Modified Gram-Schmidt Orthogonalization";
510: } else {
511: cstr = "unknown orthogonalization";
512: }
513: if (iascii) {
514: PetscViewerASCIIPrintf(viewer," restart=%D, using %s\n",gmres->max_k,cstr);
515: PetscViewerASCIIPrintf(viewer," happy breakdown tolerance %g\n",(double)gmres->haptol);
516: } else if (isstring) {
517: PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
518: }
519: return(0);
520: }
522: /*@C
523: KSPGMRESMonitorKrylov - Calls VecView() for each new direction in the GMRES accumulated Krylov space.
525: Collective on KSP
527: Input Parameters:
528: + ksp - the KSP context
529: . its - iteration number
530: . fgnorm - 2-norm of residual (or gradient)
531: - dummy - an collection of viewers created with KSPViewerCreate()
533: Options Database Keys:
534: . -ksp_gmres_kyrlov_monitor
536: Notes:
537: A new PETSCVIEWERDRAW is created for each Krylov vector so they can all be simultaneously viewed
538: Level: intermediate
540: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
542: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), KSPViewersCreate(), KSPViewersDestroy()
543: @*/
544: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
545: {
546: PetscViewers viewers = (PetscViewers)dummy;
547: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
549: Vec x;
550: PetscViewer viewer;
551: PetscBool flg;
554: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
555: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
556: if (!flg) {
557: PetscViewerSetType(viewer,PETSCVIEWERDRAW);
558: PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
559: }
560: x = VEC_VV(gmres->it+1);
561: VecView(x,viewer);
562: return(0);
563: }
565: PetscErrorCode KSPSetFromOptions_GMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
566: {
568: PetscInt restart;
569: PetscReal haptol;
570: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
571: PetscBool flg;
574: PetscOptionsHead(PetscOptionsObject,"KSP GMRES Options");
575: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
576: if (flg) { KSPGMRESSetRestart(ksp,restart); }
577: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
578: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
579: flg = PETSC_FALSE;
580: PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,NULL);
581: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
582: PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
583: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
584: PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
585: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
586: PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
587: KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
588: flg = PETSC_FALSE;
589: PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,NULL);
590: if (flg) {
591: PetscViewers viewers;
592: PetscViewersCreate(PetscObjectComm((PetscObject)ksp),&viewers);
593: KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
594: }
595: PetscOptionsTail();
596: return(0);
597: }
599: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
600: {
601: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
604: if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
605: gmres->haptol = tol;
606: return(0);
607: }
609: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
610: {
611: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
614: *max_k = gmres->max_k;
615: return(0);
616: }
618: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
619: {
620: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
624: if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
625: if (!ksp->setupstage) {
626: gmres->max_k = max_k;
627: } else if (gmres->max_k != max_k) {
628: gmres->max_k = max_k;
629: ksp->setupstage = KSP_SETUP_NEW;
630: /* free the data structures, then create them again */
631: KSPReset_GMRES(ksp);
632: }
633: return(0);
634: }
636: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
637: {
639: ((KSP_GMRES*)ksp->data)->orthog = fcn;
640: return(0);
641: }
643: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
644: {
646: *fcn = ((KSP_GMRES*)ksp->data)->orthog;
647: return(0);
648: }
650: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
651: {
652: KSP_GMRES *gmres;
655: gmres = (KSP_GMRES*)ksp->data;
656: gmres->q_preallocate = 1;
657: return(0);
658: }
660: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
661: {
662: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
665: gmres->cgstype = type;
666: return(0);
667: }
669: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
670: {
671: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
674: *type = gmres->cgstype;
675: return(0);
676: }
678: /*@
679: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
680: in the classical Gram Schmidt orthogonalization.
682: Logically Collective on KSP
684: Input Parameters:
685: + ksp - the Krylov space context
686: - type - the type of refinement
688: Options Database:
689: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
691: Level: intermediate
693: .keywords: KSP, GMRES, iterative refinement
695: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
696: KSPGMRESGetOrthogonalization()
697: @*/
698: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
699: {
705: PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
706: return(0);
707: }
709: /*@
710: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
711: in the classical Gram Schmidt orthogonalization.
713: Not Collective
715: Input Parameter:
716: . ksp - the Krylov space context
718: Output Parameter:
719: . type - the type of refinement
721: Options Database:
722: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always>
724: Level: intermediate
726: .keywords: KSP, GMRES, iterative refinement
728: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
729: KSPGMRESGetOrthogonalization()
730: @*/
731: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
732: {
737: PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
738: return(0);
739: }
742: /*@
743: KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.
745: Logically Collective on KSP
747: Input Parameters:
748: + ksp - the Krylov space context
749: - restart - integer restart value
751: Options Database:
752: . -ksp_gmres_restart <positive integer>
754: Note: The default value is 30.
756: Level: intermediate
758: .keywords: KSP, GMRES, restart, iterations
760: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
761: @*/
762: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
763: {
769: PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
770: return(0);
771: }
773: /*@
774: KSPGMRESGetRestart - Gets number of iterations at which GMRES, FGMRES and LGMRES restarts.
776: Not Collective
778: Input Parameter:
779: . ksp - the Krylov space context
781: Output Parameter:
782: . restart - integer restart value
784: Note: The default value is 30.
786: Level: intermediate
788: .keywords: KSP, GMRES, restart, iterations
790: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
791: @*/
792: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
793: {
797: PetscUseMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
798: return(0);
799: }
801: /*@
802: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
804: Logically Collective on KSP
806: Input Parameters:
807: + ksp - the Krylov space context
808: - tol - the tolerance
810: Options Database:
811: . -ksp_gmres_haptol <positive real value>
813: Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
814: a certain number of iterations. If you attempt more iterations after this point unstable
815: things can happen hence very occasionally you may need to set this value to detect this condition
817: Level: intermediate
819: .keywords: KSP, GMRES, tolerance
821: .seealso: KSPSetTolerances()
822: @*/
823: PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
824: {
829: PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
830: return(0);
831: }
833: /*MC
834: KSPGMRES - Implements the Generalized Minimal Residual method.
835: (Saad and Schultz, 1986) with restart
838: Options Database Keys:
839: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
840: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
841: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
842: vectors are allocated as needed)
843: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
844: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
845: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
846: stability of the classical Gram-Schmidt orthogonalization.
847: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
849: Level: beginner
851: Notes:
852: Left and right preconditioning are supported, but not symmetric preconditioning.
854: References:
855: . 1. - YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
856: SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.
858: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
859: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
860: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
861: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()
863: M*/
865: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
866: {
867: KSP_GMRES *gmres;
871: PetscNewLog(ksp,&gmres);
872: ksp->data = (void*)gmres;
874: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,4);
875: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
876: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,2);
877: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
878: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
880: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
881: ksp->ops->setup = KSPSetUp_GMRES;
882: ksp->ops->solve = KSPSolve_GMRES;
883: ksp->ops->reset = KSPReset_GMRES;
884: ksp->ops->destroy = KSPDestroy_GMRES;
885: ksp->ops->view = KSPView_GMRES;
886: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
887: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
888: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
889: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_ESSL)
890: ksp->ops->computeritz = KSPComputeRitz_GMRES;
891: #endif
892: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
893: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
894: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
895: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
896: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
897: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
898: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
899: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
901: gmres->haptol = 1.0e-30;
902: gmres->q_preallocate = 0;
903: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
904: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
905: gmres->nrs = 0;
906: gmres->sol_temp = 0;
907: gmres->max_k = GMRES_DEFAULT_MAXK;
908: gmres->Rsvd = 0;
909: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
910: gmres->orthogwork = 0;
911: return(0);
912: }