Actual source code: gmres.c
petsc-3.3-p7 2013-05-11
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> /*I "petscksp.h" I*/
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);
39: PetscErrorCode KSPSetUp_GMRES(KSP ksp)
40: {
41: PetscInt hh,hes,rs,cc;
43: PetscInt max_k,k;
44: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
47: max_k = gmres->max_k; /* restart size */
48: hh = (max_k + 2) * (max_k + 1);
49: hes = (max_k + 1) * (max_k + 1);
50: rs = (max_k + 2);
51: cc = (max_k + 1);
53: PetscMalloc5(hh,PetscScalar,&gmres->hh_origin,hes,PetscScalar,&gmres->hes_origin,rs,PetscScalar,&gmres->rs_origin,cc,PetscScalar,&gmres->cc_origin,cc,PetscScalar,& gmres->ss_origin);
54: PetscMemzero(gmres->hh_origin,hh*sizeof(PetscScalar));
55: PetscMemzero(gmres->hes_origin,hes*sizeof(PetscScalar));
56: PetscMemzero(gmres->rs_origin,rs*sizeof(PetscScalar));
57: PetscMemzero(gmres->cc_origin,cc*sizeof(PetscScalar));
58: PetscMemzero(gmres->ss_origin,cc*sizeof(PetscScalar));
59: PetscLogObjectMemory(ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));
61: if (ksp->calc_sings) {
62: /* Allocate workspace to hold Hessenberg matrix needed by lapack */
63: PetscMalloc((max_k + 3)*(max_k + 9)*sizeof(PetscScalar),&gmres->Rsvd);
64: PetscLogObjectMemory(ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
65: PetscMalloc(6*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
66: PetscLogObjectMemory(ksp,6*(max_k+2)*sizeof(PetscReal));
67: }
69: /* Allocate array to hold pointers to user vectors. Note that we need
70: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
71: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k + gmres->nextra_vecs;
72: PetscMalloc((gmres->vecs_allocated)*sizeof(Vec),&gmres->vecs);
73: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(Vec*),&gmres->user_work);
74: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&gmres->mwork_alloc);
75: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(sizeof(Vec*)+sizeof(PetscInt)) + gmres->vecs_allocated*sizeof(Vec));
77: if (gmres->q_preallocate) {
78: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
79: KSPGetVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,PETSC_NULL);
80: PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
81: gmres->mwork_alloc[0] = gmres->vv_allocated;
82: gmres->nwork_alloc = 1;
83: for (k=0; k<gmres->vv_allocated; k++) {
84: gmres->vecs[k] = gmres->user_work[0][k];
85: }
86: } else {
87: gmres->vv_allocated = 5;
88: KSPGetVecs(ksp,5,&gmres->user_work[0],0,PETSC_NULL);
89: PetscLogObjectParents(ksp,5,gmres->user_work[0]);
90: gmres->mwork_alloc[0] = 5;
91: gmres->nwork_alloc = 1;
92: for (k=0; k<gmres->vv_allocated; k++) {
93: gmres->vecs[k] = gmres->user_work[0][k];
94: }
95: }
96: return(0);
97: }
99: /*
100: Run gmres, possibly with restart. Return residual history if requested.
101: input parameters:
103: . gmres - structure containing parameters and work areas
105: output parameters:
106: . nres - residuals (from preconditioned system) at each step.
107: If restarting, consider passing nres+it. If null,
108: ignored
109: . itcount - number of iterations used. nres[0] to nres[itcount]
110: are defined. If null, ignored.
111:
112: Notes:
113: On entry, the value in vector VEC_VV(0) should be the initial residual
114: (this allows shortcuts where the initial preconditioned residual is 0).
115: */
118: PetscErrorCode KSPGMRESCycle(PetscInt *itcount,KSP ksp)
119: {
120: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
121: PetscReal res_norm,res,hapbnd,tt;
123: PetscInt it = 0, max_k = gmres->max_k;
124: PetscBool hapend = PETSC_FALSE;
127: VecNormalize(VEC_VV(0),&res_norm);
128: res = res_norm;
129: *GRS(0) = res_norm;
131: /* check for the convergence */
132: PetscObjectTakeAccess(ksp);
133: ksp->rnorm = res;
134: PetscObjectGrantAccess(ksp);
135: gmres->it = (it - 1);
136: KSPLogResidualHistory(ksp,res);
137: KSPMonitor(ksp,ksp->its,res);
138: if (!res) {
139: if (itcount) *itcount = 0;
140: ksp->reason = KSP_CONVERGED_ATOL;
141: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
142: return(0);
143: }
145: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
146: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
147: if (it) {
148: KSPLogResidualHistory(ksp,res);
149: KSPMonitor(ksp,ksp->its,res);
150: }
151: gmres->it = (it - 1);
152: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
153: KSPGMRESGetNewVectors(ksp,it+1);
154: }
155: KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
157: /* update hessenberg matrix and do Gram-Schmidt */
158: (*gmres->orthog)(ksp,it);
160: /* vv(i+1) . vv(i+1) */
161: VecNormalize(VEC_VV(it+1),&tt);
163: /* save the magnitude */
164: *HH(it+1,it) = tt;
165: *HES(it+1,it) = tt;
167: /* check for the happy breakdown */
168: hapbnd = PetscAbsScalar(tt / *GRS(it));
169: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
170: if (tt < hapbnd) {
171: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %14.12e tt = %14.12e\n",(double)hapbnd,(double)tt);
172: hapend = PETSC_TRUE;
173: }
174: KSPGMRESUpdateHessenberg(ksp,it,hapend,&res);
176: it++;
177: gmres->it = (it-1); /* For converged */
178: ksp->its++;
179: ksp->rnorm = res;
180: if (ksp->reason) break;
182: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
184: /* Catch error in happy breakdown and signal convergence and break from loop */
185: if (hapend) {
186: if (!ksp->reason) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down, but convergence was not indicated. Residual norm = %G",res);
187: break;
188: }
189: }
191: /* Monitor if we know that we will not return for a restart */
192: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
193: KSPLogResidualHistory(ksp,res);
194: KSPMonitor(ksp,ksp->its,res);
195: }
197: if (itcount) *itcount = it;
200: /*
201: Down here we have to solve for the "best" coefficients of the Krylov
202: columns, add the solution values together, and possibly unwind the
203: preconditioning from the solution
204: */
205: /* Form the solution (or the solution so far) */
206: KSPGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
208: return(0);
209: }
213: PetscErrorCode KSPSolve_GMRES(KSP ksp)
214: {
216: PetscInt its,itcount;
217: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
218: PetscBool guess_zero = ksp->guess_zero;
221: if (ksp->calc_sings && !gmres->Rsvd) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
223: PetscObjectTakeAccess(ksp);
224: ksp->its = 0;
225: PetscObjectGrantAccess(ksp);
227: itcount = 0;
228: ksp->reason = KSP_CONVERGED_ITERATING;
229: while (!ksp->reason) {
230: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
231: KSPGMRESCycle(&its,ksp);
232: itcount += its;
233: if (itcount >= ksp->max_it) {
234: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
235: break;
236: }
237: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
238: }
239: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
240: return(0);
241: }
245: PetscErrorCode KSPReset_GMRES(KSP ksp)
246: {
247: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
249: PetscInt i;
252: /* Free the Hessenberg matrices */
253: PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin);
255: /* free work vectors */
256: PetscFree(gmres->vecs);
257: for (i=0; i<gmres->nwork_alloc; i++) {
258: VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
259: }
260: gmres->nwork_alloc = 0;
261: PetscFree(gmres->user_work);
262: PetscFree(gmres->mwork_alloc);
263: PetscFree(gmres->nrs);
264: VecDestroy(&gmres->sol_temp);
265: PetscFree(gmres->Rsvd);
266: PetscFree(gmres->Dsvd);
267: PetscFree(gmres->orthogwork);
268: gmres->sol_temp = 0;
269: gmres->vv_allocated = 0;
270: gmres->vecs_allocated = 0;
271: gmres->sol_temp = 0;
272: return(0);
273: }
277: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
278: {
282: KSPReset_GMRES(ksp);
283: PetscFree(ksp->data);
284: /* clear composed functions */
285: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C","",PETSC_NULL);
286: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C","",PETSC_NULL);
287: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C","",PETSC_NULL);
288: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C","",PETSC_NULL);
289: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C","",PETSC_NULL);
290: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C","",PETSC_NULL);
291: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C","",PETSC_NULL);
292: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C","",PETSC_NULL);
293: return(0);
294: }
295: /*
296: KSPGMRESBuildSoln - create the solution from the starting vector and the
297: current iterates.
299: Input parameters:
300: nrs - work area of size it + 1.
301: vs - index of initial guess
302: vdest - index of result. Note that vs may == vdest (replace
303: guess with the solution).
305: This is an internal routine that knows about the GMRES internals.
306: */
309: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
310: {
311: PetscScalar tt;
313: PetscInt ii,k,j;
314: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
317: /* Solve for solution vector that minimizes the residual */
319: /* If it is < 0, no gmres steps have been performed */
320: if (it < 0) {
321: VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
322: return(0);
323: }
324: if (*HH(it,it) != 0.0) {
325: nrs[it] = *GRS(it) / *HH(it,it);
326: } else {
327: ksp->reason = KSP_DIVERGED_BREAKDOWN;
328: PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D GRS(it) = %G",it,PetscAbsScalar(*GRS(it)));
329: return(0);
330: }
331: for (ii=1; ii<=it; ii++) {
332: k = it - ii;
333: tt = *GRS(k);
334: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
335: if (*HH(k,k) == 0.0) {
336: ksp->reason = KSP_DIVERGED_BREAKDOWN;
337: PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D",k);
338: return(0);
339: }
340: nrs[k] = tt / *HH(k,k);
341: }
343: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
344: VecSet(VEC_TEMP,0.0);
345: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
347: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
348: /* add solution to previous solution */
349: if (vdest != vs) {
350: VecCopy(vs,vdest);
351: }
352: VecAXPY(vdest,1.0,VEC_TEMP);
353: return(0);
354: }
355: /*
356: Do the scalar work for the orthogonalization. Return new residual norm.
357: */
360: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
361: {
362: PetscScalar *hh,*cc,*ss,tt;
363: PetscInt j;
364: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
367: hh = HH(0,it);
368: cc = CC(0);
369: ss = SS(0);
371: /* Apply all the previously computed plane rotations to the new column
372: of the Hessenberg matrix */
373: for (j=1; j<=it; j++) {
374: tt = *hh;
375: #if defined(PETSC_USE_COMPLEX)
376: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
377: #else
378: *hh = *cc * tt + *ss * *(hh+1);
379: #endif
380: hh++;
381: *hh = *cc++ * *hh - (*ss++ * tt);
382: }
384: /*
385: compute the new plane rotation, and apply it to:
386: 1) the right-hand-side of the Hessenberg system
387: 2) the new column of the Hessenberg matrix
388: thus obtaining the updated value of the residual
389: */
390: if (!hapend) {
391: #if defined(PETSC_USE_COMPLEX)
392: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
393: #else
394: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
395: #endif
396: if (tt == 0.0) {
397: ksp->reason = KSP_DIVERGED_NULL;
398: return(0);
399: }
400: *cc = *hh / tt;
401: *ss = *(hh+1) / tt;
402: *GRS(it+1) = - (*ss * *GRS(it));
403: #if defined(PETSC_USE_COMPLEX)
404: *GRS(it) = PetscConj(*cc) * *GRS(it);
405: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
406: #else
407: *GRS(it) = *cc * *GRS(it);
408: *hh = *cc * *hh + *ss * *(hh+1);
409: #endif
410: *res = PetscAbsScalar(*GRS(it+1));
411: } else {
412: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
413: another rotation matrix (so RH doesn't change). The new residual is
414: always the new sine term times the residual from last time (GRS(it)),
415: but now the new sine rotation would be zero...so the residual should
416: be zero...so we will multiply "zero" by the last residual. This might
417: not be exactly what we want to do here -could just return "zero". */
418:
419: *res = 0.0;
420: }
421: return(0);
422: }
423: /*
424: This routine allocates more work vectors, starting from VEC_VV(it).
425: */
428: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
429: {
430: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
432: PetscInt nwork = gmres->nwork_alloc,k,nalloc;
435: nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
436: /* Adjust the number to allocate to make sure that we don't exceed the
437: number of available slots */
438: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
439: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
440: }
441: if (!nalloc) return(0);
443: gmres->vv_allocated += nalloc;
444: KSPGetVecs(ksp,nalloc,&gmres->user_work[nwork],0,PETSC_NULL);
445: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
446: gmres->mwork_alloc[nwork] = nalloc;
447: for (k=0; k<nalloc; k++) {
448: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
449: }
450: gmres->nwork_alloc++;
451: return(0);
452: }
456: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
457: {
458: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
462: if (!ptr) {
463: if (!gmres->sol_temp) {
464: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
465: PetscLogObjectParent(ksp,gmres->sol_temp);
466: }
467: ptr = gmres->sol_temp;
468: }
469: if (!gmres->nrs) {
470: /* allocate the work area */
471: PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
472: PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
473: }
475: KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
476: if (result) *result = ptr;
477: return(0);
478: }
482: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
483: {
484: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
485: const char *cstr;
487: PetscBool iascii,isstring;
490: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
491: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
492: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
493: switch (gmres->cgstype) {
494: case (KSP_GMRES_CGS_REFINE_NEVER):
495: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
496: break;
497: case (KSP_GMRES_CGS_REFINE_ALWAYS):
498: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
499: break;
500: case (KSP_GMRES_CGS_REFINE_IFNEEDED):
501: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
502: break;
503: default:
504: SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
505: }
506: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
507: cstr = "Modified Gram-Schmidt Orthogonalization";
508: } else {
509: cstr = "unknown orthogonalization";
510: }
511: if (iascii) {
512: PetscViewerASCIIPrintf(viewer," GMRES: restart=%D, using %s\n",gmres->max_k,cstr);
513: PetscViewerASCIIPrintf(viewer," GMRES: happy breakdown tolerance %G\n",gmres->haptol);
514: } else if (isstring) {
515: PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
516: } else {
517: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
518: }
519: return(0);
520: }
524: /*@C
525: KSPGMRESMonitorKrylov - Calls VecView() for each direction in the
526: GMRES accumulated Krylov space.
528: Collective on KSP
530: Input Parameters:
531: + ksp - the KSP context
532: . its - iteration number
533: . fgnorm - 2-norm of residual (or gradient)
534: - a viewers object created with PetscViewersCreate()
536: Level: intermediate
538: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
540: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
541: @*/
542: PetscErrorCode KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
543: {
544: PetscViewers viewers = (PetscViewers)dummy;
545: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
547: Vec x;
548: PetscViewer viewer;
549: PetscBool flg;
552: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
553: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
554: if (!flg) {
555: PetscViewerSetType(viewer,PETSCVIEWERDRAW);
556: PetscViewerDrawSetInfo(viewer,PETSC_NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
557: }
559: x = VEC_VV(gmres->it+1);
560: VecView(x,viewer);
562: return(0);
563: }
567: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp)
568: {
570: PetscInt restart;
571: PetscReal haptol;
572: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
573: PetscBool flg;
576: PetscOptionsHead("KSP GMRES Options");
577: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
578: if (flg) { KSPGMRESSetRestart(ksp,restart); }
579: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
580: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
581: flg = PETSC_FALSE;
582: PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,PETSC_NULL);
583: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
584: PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
585: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
586: PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
587: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
588: PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
589: KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
590: flg = PETSC_FALSE;
591: PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,PETSC_NULL);
592: if (flg) {
593: PetscViewers viewers;
594: PetscViewersCreate(((PetscObject)ksp)->comm,&viewers);
595: KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
596: }
597: PetscOptionsTail();
598: return(0);
599: }
601: extern PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
602: extern PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
605: EXTERN_C_BEGIN
608: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
609: {
610: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
613: if (tol < 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
614: gmres->haptol = tol;
615: return(0);
616: }
617: EXTERN_C_END
619: EXTERN_C_BEGIN
622: PetscErrorCode KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
623: {
624: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
627: *max_k = gmres->max_k;
628: return(0);
629: }
630: EXTERN_C_END
632: EXTERN_C_BEGIN
635: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
636: {
637: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
641: if (max_k < 1) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
642: if (!ksp->setupstage) {
643: gmres->max_k = max_k;
644: } else if (gmres->max_k != max_k) {
645: gmres->max_k = max_k;
646: ksp->setupstage = KSP_SETUP_NEW;
647: /* free the data structures, then create them again */
648: KSPReset_GMRES(ksp);
649: }
650: return(0);
651: }
652: EXTERN_C_END
654: EXTERN_C_BEGIN
657: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
658: {
660: ((KSP_GMRES *)ksp->data)->orthog = fcn;
661: return(0);
662: }
663: EXTERN_C_END
665: EXTERN_C_BEGIN
668: PetscErrorCode KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
669: {
671: *fcn = ((KSP_GMRES *)ksp->data)->orthog;
672: return(0);
673: }
674: EXTERN_C_END
676: EXTERN_C_BEGIN
679: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
680: {
681: KSP_GMRES *gmres;
684: gmres = (KSP_GMRES *)ksp->data;
685: gmres->q_preallocate = 1;
686: return(0);
687: }
688: EXTERN_C_END
690: EXTERN_C_BEGIN
693: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
694: {
695: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
698: gmres->cgstype = type;
699: return(0);
700: }
701: EXTERN_C_END
703: EXTERN_C_BEGIN
706: PetscErrorCode KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
707: {
708: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
711: *type = gmres->cgstype;
712: return(0);
713: }
714: EXTERN_C_END
718: /*@
719: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
720: in the classical Gram Schmidt orthogonalization.
722: Logically Collective on KSP
724: Input Parameters:
725: + ksp - the Krylov space context
726: - type - the type of refinement
728: Options Database:
729: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always>
731: Level: intermediate
733: .keywords: KSP, GMRES, iterative refinement
735: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
736: KSPGMRESGetOrthogonalization()
737: @*/
738: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
739: {
745: PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
746: return(0);
747: }
751: /*@
752: KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
753: in the classical Gram Schmidt orthogonalization.
755: Not Collective
757: Input Parameter:
758: . ksp - the Krylov space context
760: Output Parameter:
761: . type - the type of refinement
763: Options Database:
764: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always>
766: Level: intermediate
768: .keywords: KSP, GMRES, iterative refinement
770: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
771: KSPGMRESGetOrthogonalization()
772: @*/
773: PetscErrorCode KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
774: {
779: PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType *),(ksp,type));
780: return(0);
781: }
786: /*@
787: KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.
789: Logically Collective on KSP
791: Input Parameters:
792: + ksp - the Krylov space context
793: - restart - integer restart value
795: Options Database:
796: . -ksp_gmres_restart <positive integer>
798: Note: The default value is 30.
800: Level: intermediate
802: .keywords: KSP, GMRES, restart, iterations
804: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
805: @*/
806: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
807: {
813: PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
814: return(0);
815: }
819: /*@
820: KSPGMRESGetRestart - Gets number of iterations at which GMRES, FGMRES and LGMRES restarts.
822: Not Collective
824: Input Parameter:
825: . ksp - the Krylov space context
827: Output Parameter:
828: . restart - integer restart value
830: Note: The default value is 30.
832: Level: intermediate
834: .keywords: KSP, GMRES, restart, iterations
836: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
837: @*/
838: PetscErrorCode KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
839: {
843: PetscTryMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
844: return(0);
845: }
849: /*@
850: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
852: Logically Collective on KSP
854: Input Parameters:
855: + ksp - the Krylov space context
856: - tol - the tolerance
858: Options Database:
859: . -ksp_gmres_haptol <positive real value>
861: Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
862: a certain number of iterations. If you attempt more iterations after this point unstable
863: things can happen hence very occasionally you may need to set this value to detect this condition
865: Level: intermediate
867: .keywords: KSP, GMRES, tolerance
869: .seealso: KSPSetTolerances()
870: @*/
871: PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
872: {
877: PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
878: return(0);
879: }
881: /*MC
882: KSPGMRES - Implements the Generalized Minimal Residual method.
883: (Saad and Schultz, 1986) with restart
886: Options Database Keys:
887: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
888: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
889: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
890: vectors are allocated as needed)
891: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
892: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
893: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
894: stability of the classical Gram-Schmidt orthogonalization.
895: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
897: Level: beginner
899: Notes: Left and right preconditioning are supported, but not symmetric preconditioning.
901: References:
902: GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS. YOUCEF SAAD AND MARTIN H. SCHULTZ,
903: SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986, pp. 856--869.
905: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
906: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
907: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
908: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()
910: M*/
912: EXTERN_C_BEGIN
915: PetscErrorCode KSPCreate_GMRES(KSP ksp)
916: {
917: KSP_GMRES *gmres;
921: PetscNewLog(ksp,KSP_GMRES,&gmres);
922: ksp->data = (void*)gmres;
924: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
925: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
927: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
928: ksp->ops->setup = KSPSetUp_GMRES;
929: ksp->ops->solve = KSPSolve_GMRES;
930: ksp->ops->reset = KSPReset_GMRES;
931: ksp->ops->destroy = KSPDestroy_GMRES;
932: ksp->ops->view = KSPView_GMRES;
933: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
934: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
935: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
937: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
938: "KSPGMRESSetPreAllocateVectors_GMRES",
939: KSPGMRESSetPreAllocateVectors_GMRES);
940: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
941: "KSPGMRESSetOrthogonalization_GMRES",
942: KSPGMRESSetOrthogonalization_GMRES);
943: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
944: "KSPGMRESGetOrthogonalization_GMRES",
945: KSPGMRESGetOrthogonalization_GMRES);
946: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
947: "KSPGMRESSetRestart_GMRES",
948: KSPGMRESSetRestart_GMRES);
949: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
950: "KSPGMRESGetRestart_GMRES",
951: KSPGMRESGetRestart_GMRES);
952: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
953: "KSPGMRESSetHapTol_GMRES",
954: KSPGMRESSetHapTol_GMRES);
955: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
956: "KSPGMRESSetCGSRefinementType_GMRES",
957: KSPGMRESSetCGSRefinementType_GMRES);
958: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
959: "KSPGMRESGetCGSRefinementType_GMRES",
960: KSPGMRESGetCGSRefinementType_GMRES);
962: gmres->haptol = 1.0e-30;
963: gmres->q_preallocate = 0;
964: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
965: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
966: gmres->nrs = 0;
967: gmres->sol_temp = 0;
968: gmres->max_k = GMRES_DEFAULT_MAXK;
969: gmres->Rsvd = 0;
970: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
971: gmres->orthogwork = 0;
972: return(0);
973: }
974: EXTERN_C_END