Actual source code: gmres.c
1: /*
2: This file implements GMRES (a Generalized Minimal Residual) method.
3: Reference: Saad and Schultz, 1986.
6: Some comments on left vs. right preconditioning, and restarts.
7: Left and right preconditioning.
8: If right preconditioning is chosen, then the problem being solved
9: by gmres is actually
10: My = AB^-1 y = f
11: so the initial residual is
12: r = f - Mx
13: Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
14: residual is
15: r = f - A x
16: The final solution is then
17: x = B^-1 y
19: If left preconditioning is chosen, then the problem being solved is
20: My = B^-1 A x = B^-1 f,
21: and the initial residual is
22: r = B^-1(f - Ax)
24: Restarts: Restarts are basically solves with x0 not equal to zero.
25: Note that we can eliminate an extra application of B^-1 between
26: restarts as long as we don't require that the solution at the end
27: of an unsuccessful gmres iteration always be the solution x.
28: */
30: #include src/ksp/ksp/impls/gmres/gmresp.h
31: #define GMRES_DELTA_DIRECTIONS 10
32: #define GMRES_DEFAULT_MAXK 30
33: static PetscErrorCode GMRESGetNewVectors(KSP,PetscInt);
34: static PetscErrorCode GMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal*);
35: static PetscErrorCode BuildGmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
39: PetscErrorCode KSPSetUp_GMRES(KSP ksp)
40: {
41: PetscInt size,hh,hes,rs,cc;
43: PetscInt max_k,k;
44: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
45: Vec vec;
46: Mat pmat;
49: if (ksp->pc_side == PC_SYMMETRIC) {
50: SETERRQ(2,"no symmetric preconditioning for KSPGMRES");
51: }
53: max_k = gmres->max_k;
54: hh = (max_k + 2) * (max_k + 1);
55: hes = (max_k + 1) * (max_k + 1);
56: rs = (max_k + 2);
57: cc = (max_k + 1);
58: size = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);
60: PetscMalloc(size,&gmres->hh_origin);
61: PetscMemzero(gmres->hh_origin,size);
62: PetscLogObjectMemory(ksp,size);
63: gmres->hes_origin = gmres->hh_origin + hh;
64: gmres->rs_origin = gmres->hes_origin + hes;
65: gmres->cc_origin = gmres->rs_origin + rs;
66: gmres->ss_origin = gmres->cc_origin + cc;
68: if (ksp->calc_sings) {
69: /* Allocate workspace to hold Hessenberg matrix needed by lapack */
70: size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
71: PetscMalloc(size,&gmres->Rsvd);
72: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
73: PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
74: }
76: /* Allocate array to hold pointers to user vectors. Note that we need
77: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
78: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&gmres->vecs);
79: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
80: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&gmres->user_work);
81: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&gmres->mwork_alloc);
82: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));
84: PCGetOperators(ksp->pc,0,&pmat,0);
85: if (!pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
86: MatGetVecs(pmat,&vec,0);
87: if (gmres->q_preallocate) {
88: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
89: KSPGetVecs(ksp,gmres->vv_allocated,&gmres->user_work[0]);
90: PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
91: gmres->mwork_alloc[0] = gmres->vv_allocated;
92: gmres->nwork_alloc = 1;
93: for (k=0; k<gmres->vv_allocated; k++) {
94: gmres->vecs[k] = gmres->user_work[0][k];
95: }
96: } else {
97: gmres->vv_allocated = 5;
98: KSPGetVecs(ksp,5,&gmres->user_work[0]);
99: PetscLogObjectParents(ksp,5,gmres->user_work[0]);
100: gmres->mwork_alloc[0] = 5;
101: gmres->nwork_alloc = 1;
102: for (k=0; k<gmres->vv_allocated; k++) {
103: gmres->vecs[k] = gmres->user_work[0][k];
104: }
105: }
106: VecDestroy(vec);
107: return(0);
108: }
110: /*
111: Run gmres, possibly with restart. Return residual history if requested.
112: input parameters:
114: . gmres - structure containing parameters and work areas
116: output parameters:
117: . nres - residuals (from preconditioned system) at each step.
118: If restarting, consider passing nres+it. If null,
119: ignored
120: . itcount - number of iterations used. nres[0] to nres[itcount]
121: are defined. If null, ignored.
122:
123: Notes:
124: On entry, the value in vector VEC_VV(0) should be the initial residual
125: (this allows shortcuts where the initial preconditioned residual is 0).
126: */
129: PetscErrorCode GMREScycle(PetscInt *itcount,KSP ksp)
130: {
131: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
132: PetscReal res_norm,res,hapbnd,tt;
134: PetscInt it = 0, max_k = gmres->max_k;
135: PetscTruth hapend = PETSC_FALSE;
138: VecNormalize(VEC_VV(0),&res_norm);
139: res = res_norm;
140: *GRS(0) = res_norm;
142: /* check for the convergence */
143: PetscObjectTakeAccess(ksp);
144: ksp->rnorm = res;
145: PetscObjectGrantAccess(ksp);
146: gmres->it = (it - 1);
147: KSPLogResidualHistory(ksp,res);
148: if (!res) {
149: if (itcount) *itcount = 0;
150: ksp->reason = KSP_CONVERGED_ATOL;
151: PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entry\n");
152: return(0);
153: }
155: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
156: while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
157: KSPLogResidualHistory(ksp,res);
158: gmres->it = (it - 1);
159: KSPMonitor(ksp,ksp->its,res);
160: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
161: GMRESGetNewVectors(ksp,it+1);
162: }
163: KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
165: /* update hessenberg matrix and do Gram-Schmidt */
166: (*gmres->orthog)(ksp,it);
168: /* vv(i+1) . vv(i+1) */
169: VecNormalize(VEC_VV(it+1),&tt);
170: /* save the magnitude */
171: *HH(it+1,it) = tt;
172: *HES(it+1,it) = tt;
174: /* check for the happy breakdown */
175: hapbnd = PetscAbsScalar(tt / *GRS(it));
176: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
177: if (tt < hapbnd) {
178: PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",hapbnd,tt);
179: hapend = PETSC_TRUE;
180: }
181: GMRESUpdateHessenberg(ksp,it,hapend,&res);
182: if (ksp->reason) break;
184: it++;
185: gmres->it = (it-1); /* For converged */
186: PetscObjectTakeAccess(ksp);
187: ksp->its++;
188: ksp->rnorm = res;
189: PetscObjectGrantAccess(ksp);
191: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
193: /* Catch error in happy breakdown and signal convergence and break from loop */
194: if (hapend) {
195: if (!ksp->reason) {
196: SETERRQ1(0,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",res);
197: }
198: break;
199: }
200: }
202: /* Monitor if we know that we will not return for a restart */
203: if (ksp->reason || ksp->its >= ksp->max_it) {
204: KSPLogResidualHistory(ksp,res);
205: KSPMonitor(ksp,ksp->its,res);
206: }
208: if (itcount) *itcount = it;
211: /*
212: Down here we have to solve for the "best" coefficients of the Krylov
213: columns, add the solution values together, and possibly unwind the
214: preconditioning from the solution
215: */
216: /* Form the solution (or the solution so far) */
217: BuildGmresSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
219: return(0);
220: }
224: PetscErrorCode KSPSolve_GMRES(KSP ksp)
225: {
227: PetscInt its,itcount;
228: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
229: PetscTruth guess_zero = ksp->guess_zero;
232: if (ksp->calc_sings && !gmres->Rsvd) {
233: SETERRQ(PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
234: }
236: PetscObjectTakeAccess(ksp);
237: ksp->its = 0;
238: PetscObjectGrantAccess(ksp);
240: itcount = 0;
241: ksp->reason = KSP_CONVERGED_ITERATING;
242: while (!ksp->reason) {
243: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
244: GMREScycle(&its,ksp);
245: itcount += its;
246: if (itcount >= ksp->max_it) {
247: ksp->reason = KSP_DIVERGED_ITS;
248: break;
249: }
250: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
251: }
252: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
253: return(0);
254: }
258: PetscErrorCode KSPDestroy_GMRES_Internal(KSP ksp)
259: {
260: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
262: PetscInt i;
265: /* Free the Hessenberg matrix */
266: if (gmres->hh_origin) {
267: PetscFree(gmres->hh_origin);
268: gmres->hh_origin = 0;
269: }
271: /* Free the pointer to user variables */
272: if (gmres->vecs) {
273: PetscFree(gmres->vecs);
274: gmres->vecs = 0;
275: }
277: /* free work vectors */
278: for (i=0; i<gmres->nwork_alloc; i++) {
279: VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
280: }
281: if (gmres->user_work) {
282: PetscFree(gmres->user_work);
283: gmres->user_work = 0;
284: }
285: if (gmres->mwork_alloc) {
286: PetscFree(gmres->mwork_alloc);
287: gmres->mwork_alloc = 0;
288: }
289: if (gmres->nrs) {
290: PetscFree(gmres->nrs);
291: gmres->nrs = 0;
292: }
293: if (gmres->sol_temp) {
294: VecDestroy(gmres->sol_temp);
295: gmres->sol_temp = 0;
296: }
297: if (gmres->Rsvd) {
298: PetscFree(gmres->Rsvd);
299: gmres->Rsvd = 0;
300: }
301: if (gmres->Dsvd) {
302: PetscFree(gmres->Dsvd);
303: gmres->Dsvd = 0;
304: }
306: gmres->nwork_alloc = 0;
307: gmres->vv_allocated = 0;
308: gmres->vecs_allocated = 0;
309: gmres->nrs = 0;
310: gmres->sol_temp = 0;
311: gmres->Rsvd = 0;
312: return(0);
313: }
317: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
318: {
319: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
323: KSPDestroy_GMRES_Internal(ksp);
324: PetscFree(gmres);
325: return(0);
326: }
327: /*
328: BuildGmresSoln - create the solution from the starting vector and the
329: current iterates.
331: Input parameters:
332: nrs - work area of size it + 1.
333: vs - index of initial guess
334: vdest - index of result. Note that vs may == vdest (replace
335: guess with the solution).
337: This is an internal routine that knows about the GMRES internals.
338: */
341: static PetscErrorCode BuildGmresSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
342: {
343: PetscScalar tt,zero = 0.0,one = 1.0;
345: PetscInt ii,k,j;
346: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
349: /* Solve for solution vector that minimizes the residual */
351: /* If it is < 0, no gmres steps have been performed */
352: if (it < 0) {
353: if (vdest != vs) {
354: VecCopy(vs,vdest);
355: }
356: return(0);
357: }
358: if (*HH(it,it) == 0.0) SETERRQ2(PETSC_ERR_CONV_FAILED,"HH(it,it) is identically zero; it = %D GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
359: if (*HH(it,it) != 0.0) {
360: nrs[it] = *GRS(it) / *HH(it,it);
361: } else {
362: nrs[it] = 0.0;
363: }
364: for (ii=1; ii<=it; ii++) {
365: k = it - ii;
366: tt = *GRS(k);
367: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
368: nrs[k] = tt / *HH(k,k);
369: }
371: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
372: VecSet(&zero,VEC_TEMP);
373: VecMAXPY(it+1,nrs,VEC_TEMP,&VEC_VV(0));
375: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
376: /* add solution to previous solution */
377: if (vdest != vs) {
378: VecCopy(vs,vdest);
379: }
380: VecAXPY(&one,VEC_TEMP,vdest);
381: return(0);
382: }
383: /*
384: Do the scalar work for the orthogonalization. Return new residual.
385: */
388: static PetscErrorCode GMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
389: {
390: PetscScalar *hh,*cc,*ss,tt;
391: PetscInt j;
392: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
395: hh = HH(0,it);
396: cc = CC(0);
397: ss = SS(0);
399: /* Apply all the previously computed plane rotations to the new column
400: of the Hessenberg matrix */
401: for (j=1; j<=it; j++) {
402: tt = *hh;
403: #if defined(PETSC_USE_COMPLEX)
404: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
405: #else
406: *hh = *cc * tt + *ss * *(hh+1);
407: #endif
408: hh++;
409: *hh = *cc++ * *hh - (*ss++ * tt);
410: }
412: /*
413: compute the new plane rotation, and apply it to:
414: 1) the right-hand-side of the Hessenberg system
415: 2) the new column of the Hessenberg matrix
416: thus obtaining the updated value of the residual
417: */
418: if (!hapend) {
419: #if defined(PETSC_USE_COMPLEX)
420: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
421: #else
422: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
423: #endif
424: if (tt == 0.0) {
425: ksp->reason = KSP_DIVERGED_NULL;
426: return(0);
427: }
428: *cc = *hh / tt;
429: *ss = *(hh+1) / tt;
430: *GRS(it+1) = - (*ss * *GRS(it));
431: #if defined(PETSC_USE_COMPLEX)
432: *GRS(it) = PetscConj(*cc) * *GRS(it);
433: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
434: #else
435: *GRS(it) = *cc * *GRS(it);
436: *hh = *cc * *hh + *ss * *(hh+1);
437: #endif
438: *res = PetscAbsScalar(*GRS(it+1));
439: } else {
440: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
441: another rotation matrix (so RH doesn't change). The new residual is
442: always the new sine term times the residual from last time (GRS(it)),
443: but now the new sine rotation would be zero...so the residual should
444: be zero...so we will multiply "zero" by the last residual. This might
445: not be exactly what we want to do here -could just return "zero". */
446:
447: *res = 0.0;
448: }
449: return(0);
450: }
451: /*
452: This routine allocates more work vectors, starting from VEC_VV(it).
453: */
456: static PetscErrorCode GMRESGetNewVectors(KSP ksp,PetscInt it)
457: {
458: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
460: PetscInt nwork = gmres->nwork_alloc,k,nalloc;
463: nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
464: /* Adjust the number to allocate to make sure that we don't exceed the
465: number of available slots */
466: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
467: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
468: }
469: if (!nalloc) return(0);
471: gmres->vv_allocated += nalloc;
472: KSPGetVecs(ksp,nalloc,&gmres->user_work[nwork]);
473: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
474: gmres->mwork_alloc[nwork] = nalloc;
475: for (k=0; k<nalloc; k++) {
476: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
477: }
478: gmres->nwork_alloc++;
479: return(0);
480: }
484: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
485: {
486: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
490: if (!ptr) {
491: if (!gmres->sol_temp) {
492: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
493: PetscLogObjectParent(ksp,gmres->sol_temp);
494: }
495: ptr = gmres->sol_temp;
496: }
497: if (!gmres->nrs) {
498: /* allocate the work area */
499: PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
500: PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
501: }
503: BuildGmresSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
504: *result = ptr;
505: return(0);
506: }
510: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
511: {
512: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
513: const char *cstr;
515: PetscTruth iascii,isstring;
518: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
519: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
520: if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
521: if (gmres->cgstype == KSP_GMRES_CGS_REFINE_NEVER) {
522: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
523: } else if (gmres->cgstype == KSP_GMRES_CGS_REFINE_ALWAYS) {
524: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
525: } else {
526: cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
527: }
528: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
529: cstr = "Modified Gram-Schmidt Orthogonalization";
530: } else {
531: cstr = "unknown orthogonalization";
532: }
533: if (iascii) {
534: PetscViewerASCIIPrintf(viewer," GMRES: restart=%D, using %s\n",gmres->max_k,cstr);
535: PetscViewerASCIIPrintf(viewer," GMRES: happy breakdown tolerance %g\n",gmres->haptol);
536: } else if (isstring) {
537: PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
538: } else {
539: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
540: }
541: return(0);
542: }
546: /*@C
547: KSPGMRESKrylovMonitor - Calls VecView() for each direction in the
548: GMRES accumulated Krylov space.
550: Collective on KSP
552: Input Parameters:
553: + ksp - the KSP context
554: . its - iteration number
555: . fgnorm - 2-norm of residual (or gradient)
556: - a viewers object created with PetscViewersCreate()
558: Level: intermediate
560: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
562: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
563: @*/
564: PetscErrorCode KSPGMRESKrylovMonitor(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
565: {
566: PetscViewers viewers = (PetscViewers)dummy;
567: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
569: Vec x;
570: PetscViewer viewer;
573: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
574: PetscViewerSetType(viewer,PETSC_VIEWER_DRAW);
576: x = VEC_VV(gmres->it+1);
577: VecView(x,viewer);
579: return(0);
580: }
584: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp)
585: {
587: PetscInt restart,indx;
588: PetscReal haptol;
589: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
590: PetscTruth flg;
591: const char *types[] = {"never","ifneeded","always"};
594: PetscOptionsHead("KSP GMRES Options");
595: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
596: if (flg) { KSPGMRESSetRestart(ksp,restart); }
597: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
598: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
599: PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
600: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
601: PetscOptionsLogicalGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
602: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
603: PetscOptionsLogicalGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
604: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
605: PetscOptionsEList("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType()",types,3,types[(PetscInt)gmres->cgstype],&indx,&flg);
606: if (flg) {
607: KSPGMRESSetCGSRefinementType(ksp,(KSPGMRESCGSRefinementType)indx);
608: }
610: PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
611: if (flg) {
612: PetscViewers viewers;
613: PetscViewersCreate(ksp->comm,&viewers);
614: KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(PetscErrorCode (*)(void*))PetscViewersDestroy);
615: }
616: PetscOptionsTail();
617: return(0);
618: }
620: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
621: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
627: PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
628: {
629: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
632: if (tol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
633: gmres->haptol = tol;
634: return(0);
635: }
641: PetscErrorCode KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
642: {
643: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
647: if (max_k < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
648: if (!ksp->setupcalled) {
649: gmres->max_k = max_k;
650: } else if (gmres->max_k != max_k) {
651: gmres->max_k = max_k;
652: ksp->setupcalled = 0;
653: /* free the data structures, then create them again */
654: KSPDestroy_GMRES_Internal(ksp);
655: }
657: return(0);
658: }
665: PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
666: {
669: ((KSP_GMRES *)ksp->data)->orthog = fcn;
670: return(0);
671: }
677: PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
678: {
679: KSP_GMRES *gmres;
682: gmres = (KSP_GMRES *)ksp->data;
683: gmres->q_preallocate = 1;
684: return(0);
685: }
691: PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
692: {
693: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
696: gmres->cgstype = type;
697: return(0);
698: }
703: /*@
704: KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
705: in the classical Gram Schmidt orthogonalization.
706: of the preconditioned problem.
708: Collective on KSP
710: Input Parameters:
711: + ksp - the Krylov space context
712: - type - the type of refinement
714: Options Database:
715: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always>
717: Level: intermediate
719: .keywords: KSP, GMRES, iterative refinement
721: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization()
722: @*/
723: PetscErrorCode KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
724: {
725: PetscErrorCode ierr,(*f)(KSP,KSPGMRESCGSRefinementType);
729: PetscObjectQueryFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",(void (**)(void))&f);
730: if (f) {
731: (*f)(ksp,type);
732: }
733: return(0);
734: }
738: /*@C
739: KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.
741: Collective on KSP
743: Input Parameters:
744: + ksp - the Krylov space context
745: - restart - integer restart value
747: Options Database:
748: . -ksp_gmres_restart <positive integer>
750: Note: The default value is 30.
752: Level: intermediate
754: .keywords: KSP, GMRES, restart, iterations
756: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors()
757: @*/
758: PetscErrorCode KSPGMRESSetRestart(KSP ksp, PetscInt restart)
759: {
763: PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,PetscInt),((ksp),(restart)));
764: return(0);
765: }
769: /*@
770: KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.
772: Collective on KSP
774: Input Parameters:
775: + ksp - the Krylov space context
776: - tol - the tolerance
778: Options Database:
779: . -ksp_gmres_haptol <positive real value>
781: Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
782: a certain number of iterations. If you attempt more iterations after this point unstable
783: things can happen hence very occasionally you may need to set this value to detect this condition
785: Level: intermediate
787: .keywords: KSP, GMRES, tolerance
789: .seealso: KSPSetTolerances()
790: @*/
791: PetscErrorCode KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
792: {
796: PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)));
797: return(0);
798: }
800: /*MC
801: KSPGMRES - Implements the Generalized Minimal Residual method.
802: (Saad and Schultz, 1986) with restart
805: Options Database Keys:
806: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
807: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
808: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
809: vectors are allocated as needed)
810: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
811: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
812: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
813: stability of the classical Gram-Schmidt orthogonalization.
814: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
816: Level: beginner
819: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
820: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
821: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
822: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor()
824: M*/
829: PetscErrorCode KSPCreate_GMRES(KSP ksp)
830: {
831: KSP_GMRES *gmres;
835: PetscNew(KSP_GMRES,&gmres);
836: PetscMemzero(gmres,sizeof(KSP_GMRES));
837: PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
838: ksp->data = (void*)gmres;
839: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
841: ksp->ops->setup = KSPSetUp_GMRES;
842: ksp->ops->solve = KSPSolve_GMRES;
843: ksp->ops->destroy = KSPDestroy_GMRES;
844: ksp->ops->view = KSPView_GMRES;
845: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
846: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
847: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
849: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
850: "KSPGMRESSetPreAllocateVectors_GMRES",
851: KSPGMRESSetPreAllocateVectors_GMRES);
852: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
853: "KSPGMRESSetOrthogonalization_GMRES",
854: KSPGMRESSetOrthogonalization_GMRES);
855: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
856: "KSPGMRESSetRestart_GMRES",
857: KSPGMRESSetRestart_GMRES);
858: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
859: "KSPGMRESSetHapTol_GMRES",
860: KSPGMRESSetHapTol_GMRES);
861: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
862: "KSPGMRESSetCGSRefinementType_GMRES",
863: KSPGMRESSetCGSRefinementType_GMRES);
865: gmres->haptol = 1.0e-30;
866: gmres->q_preallocate = 0;
867: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
868: gmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
869: gmres->nrs = 0;
870: gmres->sol_temp = 0;
871: gmres->max_k = GMRES_DEFAULT_MAXK;
872: gmres->Rsvd = 0;
873: gmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
874: return(0);
875: }