Actual source code: gmres.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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:   PetscCalloc5(hh,&gmres->hh_origin,hes,&gmres->hes_origin,rs,&gmres->rs_origin,cc,&gmres->cc_origin,cc,&gmres->ss_origin);
 54:   PetscLogObjectMemory((PetscObject)ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));

 56:   if (ksp->calc_sings) {
 57:     /* Allocate workspace to hold Hessenberg matrix needed by lapack */
 58:     PetscMalloc1((max_k + 3)*(max_k + 9),&gmres->Rsvd);
 59:     PetscLogObjectMemory((PetscObject)ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
 60:     PetscMalloc1(6*(max_k+2),&gmres->Dsvd);
 61:     PetscLogObjectMemory((PetscObject)ksp,6*(max_k+2)*sizeof(PetscReal));
 62:   }

 64:   /* Allocate array to hold pointers to user vectors.  Note that we need
 65:    4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
 66:   gmres->vecs_allocated = VEC_OFFSET + 2 + max_k + gmres->nextra_vecs;

 68:   PetscMalloc1((gmres->vecs_allocated),&gmres->vecs);
 69:   PetscMalloc1((VEC_OFFSET+2+max_k),&gmres->user_work);
 70:   PetscMalloc1((VEC_OFFSET+2+max_k),&gmres->mwork_alloc);
 71:   PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+2+max_k)*(sizeof(Vec*)+sizeof(PetscInt)) + gmres->vecs_allocated*sizeof(Vec));

 73:   if (gmres->q_preallocate) {
 74:     gmres->vv_allocated = VEC_OFFSET + 2 + max_k;

 76:     KSPGetVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,NULL);
 77:     PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);

 79:     gmres->mwork_alloc[0] = gmres->vv_allocated;
 80:     gmres->nwork_alloc    = 1;
 81:     for (k=0; k<gmres->vv_allocated; k++) {
 82:       gmres->vecs[k] = gmres->user_work[0][k];
 83:     }
 84:   } else {
 85:     gmres->vv_allocated = 5;

 87:     KSPGetVecs(ksp,5,&gmres->user_work[0],0,NULL);
 88:     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.

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:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
133:   ksp->rnorm = res;
134:   PetscObjectSAWsGrantAccess((PetscObject)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) {
187:         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);
188:         else {
189:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
190:           break;
191:         }
192:       }
193:     }
194:   }

196:   /* Monitor if we know that we will not return for a restart */
197:   if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
198:     KSPLogResidualHistory(ksp,res);
199:     KSPMonitor(ksp,ksp->its,res);
200:   }

202:   if (itcount) *itcount = it;


205:   /*
206:     Down here we have to solve for the "best" coefficients of the Krylov
207:     columns, add the solution values together, and possibly unwind the
208:     preconditioning from the solution
209:    */
210:   /* Form the solution (or the solution so far) */
211:   KSPGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-1);
212:   return(0);
213: }

217: PetscErrorCode KSPSolve_GMRES(KSP ksp)
218: {
220:   PetscInt       its,itcount;
221:   KSP_GMRES      *gmres     = (KSP_GMRES*)ksp->data;
222:   PetscBool      guess_zero = ksp->guess_zero;

225:   if (ksp->calc_sings && !gmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");

227:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
228:   ksp->its = 0;
229:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

231:   itcount     = 0;
232:   ksp->reason = KSP_CONVERGED_ITERATING;
233:   while (!ksp->reason) {
234:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
235:     KSPGMRESCycle(&its,ksp);
236:     itcount += its;
237:     if (itcount >= ksp->max_it) {
238:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
239:       break;
240:     }
241:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
242:   }
243:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
244:   return(0);
245: }

249: PetscErrorCode KSPReset_GMRES(KSP ksp)
250: {
251:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
253:   PetscInt       i;

256:   /* Free the Hessenberg matrices */
257:   PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin);

259:   /* free work vectors */
260:   PetscFree(gmres->vecs);
261:   for (i=0; i<gmres->nwork_alloc; i++) {
262:     VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
263:   }
264:   gmres->nwork_alloc = 0;

266:   PetscFree(gmres->user_work);
267:   PetscFree(gmres->mwork_alloc);
268:   PetscFree(gmres->nrs);
269:   VecDestroy(&gmres->sol_temp);
270:   PetscFree(gmres->Rsvd);
271:   PetscFree(gmres->Dsvd);
272:   PetscFree(gmres->orthogwork);

274:   gmres->sol_temp       = 0;
275:   gmres->vv_allocated   = 0;
276:   gmres->vecs_allocated = 0;
277:   gmres->sol_temp       = 0;
278:   return(0);
279: }

283: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
284: {

288:   KSPReset_GMRES(ksp);
289:   PetscFree(ksp->data);
290:   /* clear composed functions */
291:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",NULL);
292:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",NULL);
293:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",NULL);
294:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",NULL);
295:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",NULL);
296:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",NULL);
297:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",NULL);
298:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",NULL);
299:   return(0);
300: }
301: /*
302:     KSPGMRESBuildSoln - create the solution from the starting vector and the
303:     current iterates.

305:     Input parameters:
306:         nrs - work area of size it + 1.
307:         vs  - index of initial guess
308:         vdest - index of result.  Note that vs may == vdest (replace
309:                 guess with the solution).

311:      This is an internal routine that knows about the GMRES internals.
312:  */
315: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
316: {
317:   PetscScalar    tt;
319:   PetscInt       ii,k,j;
320:   KSP_GMRES      *gmres = (KSP_GMRES*)(ksp->data);

323:   /* Solve for solution vector that minimizes the residual */

325:   /* If it is < 0, no gmres steps have been performed */
326:   if (it < 0) {
327:     VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
328:     return(0);
329:   }
330:   if (*HH(it,it) != 0.0) {
331:     nrs[it] = *GRS(it) / *HH(it,it);
332:   } else {
333:     ksp->reason = KSP_DIVERGED_BREAKDOWN;

335:     PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D GRS(it) = %g",it,(double)PetscAbsScalar(*GRS(it)));
336:     return(0);
337:   }
338:   for (ii=1; ii<=it; ii++) {
339:     k  = it - ii;
340:     tt = *GRS(k);
341:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
342:     if (*HH(k,k) == 0.0) {
343:       ksp->reason = KSP_DIVERGED_BREAKDOWN;

345:       PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D",k);
346:       return(0);
347:     }
348:     nrs[k] = tt / *HH(k,k);
349:   }

351:   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
352:   VecSet(VEC_TEMP,0.0);
353:   VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));

355:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
356:   /* add solution to previous solution */
357:   if (vdest != vs) {
358:     VecCopy(vs,vdest);
359:   }
360:   VecAXPY(vdest,1.0,VEC_TEMP);
361:   return(0);
362: }
363: /*
364:    Do the scalar work for the orthogonalization.  Return new residual norm.
365:  */
368: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
369: {
370:   PetscScalar *hh,*cc,*ss,tt;
371:   PetscInt    j;
372:   KSP_GMRES   *gmres = (KSP_GMRES*)(ksp->data);

375:   hh = HH(0,it);
376:   cc = CC(0);
377:   ss = SS(0);

379:   /* Apply all the previously computed plane rotations to the new column
380:      of the Hessenberg matrix */
381:   for (j=1; j<=it; j++) {
382:     tt  = *hh;
383:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
384:     hh++;
385:     *hh = *cc++ * *hh - (*ss++ * tt);
386:   }

388:   /*
389:     compute the new plane rotation, and apply it to:
390:      1) the right-hand-side of the Hessenberg system
391:      2) the new column of the Hessenberg matrix
392:     thus obtaining the updated value of the residual
393:   */
394:   if (!hapend) {
395:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
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:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
404:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
405:     *res       = PetscAbsScalar(*GRS(it+1));
406:   } else {
407:     /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
408:             another rotation matrix (so RH doesn't change).  The new residual is
409:             always the new sine term times the residual from last time (GRS(it)),
410:             but now the new sine rotation would be zero...so the residual should
411:             be zero...so we will multiply "zero" by the last residual.  This might
412:             not be exactly what we want to do here -could just return "zero". */

414:     *res = 0.0;
415:   }
416:   return(0);
417: }
418: /*
419:    This routine allocates more work vectors, starting from VEC_VV(it).
420:  */
423: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
424: {
425:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
427:   PetscInt       nwork = gmres->nwork_alloc,k,nalloc;

430:   nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
431:   /* Adjust the number to allocate to make sure that we don't exceed the
432:     number of available slots */
433:   if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) {
434:     nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
435:   }
436:   if (!nalloc) return(0);

438:   gmres->vv_allocated += nalloc;

440:   KSPGetVecs(ksp,nalloc,&gmres->user_work[nwork],0,NULL);
441:   PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);

443:   gmres->mwork_alloc[nwork] = nalloc;
444:   for (k=0; k<nalloc; k++) {
445:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
446:   }
447:   gmres->nwork_alloc++;
448:   return(0);
449: }

453: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
454: {
455:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

459:   if (!ptr) {
460:     if (!gmres->sol_temp) {
461:       VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
462:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)gmres->sol_temp);
463:     }
464:     ptr = gmres->sol_temp;
465:   }
466:   if (!gmres->nrs) {
467:     /* allocate the work area */
468:     PetscMalloc1(gmres->max_k,&gmres->nrs);
469:     PetscLogObjectMemory((PetscObject)ksp,gmres->max_k*sizeof(PetscScalar));
470:   }

472:   KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
473:   if (result) *result = ptr;
474:   return(0);
475: }

479: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
480: {
481:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
482:   const char     *cstr;
484:   PetscBool      iascii,isstring;

487:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
488:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
489:   if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
490:     switch (gmres->cgstype) {
491:     case (KSP_GMRES_CGS_REFINE_NEVER):
492:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
493:       break;
494:     case (KSP_GMRES_CGS_REFINE_ALWAYS):
495:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
496:       break;
497:     case (KSP_GMRES_CGS_REFINE_IFNEEDED):
498:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
499:       break;
500:     default:
501:       SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
502:     }
503:   } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
504:     cstr = "Modified Gram-Schmidt Orthogonalization";
505:   } else {
506:     cstr = "unknown orthogonalization";
507:   }
508:   if (iascii) {
509:     PetscViewerASCIIPrintf(viewer,"  GMRES: restart=%D, using %s\n",gmres->max_k,cstr);
510:     PetscViewerASCIIPrintf(viewer,"  GMRES: happy breakdown tolerance %g\n",(double)gmres->haptol);
511:   } else if (isstring) {
512:     PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
513:   }
514:   return(0);
515: }

519: /*@C
520:    KSPGMRESMonitorKrylov - Calls VecView() for each direction in the
521:    GMRES accumulated Krylov space.

523:    Collective on KSP

525:    Input Parameters:
526: +  ksp - the KSP context
527: .  its - iteration number
528: .  fgnorm - 2-norm of residual (or gradient)
529: -  a viewers object created with PetscViewersCreate()

531:    Level: intermediate

533: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space

535: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
536: @*/
537: PetscErrorCode  KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
538: {
539:   PetscViewers   viewers = (PetscViewers)dummy;
540:   KSP_GMRES      *gmres  = (KSP_GMRES*)ksp->data;
542:   Vec            x;
543:   PetscViewer    viewer;
544:   PetscBool      flg;

547:   PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
548:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
549:   if (!flg) {
550:     PetscViewerSetType(viewer,PETSCVIEWERDRAW);
551:     PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
552:   }

554:   x    = VEC_VV(gmres->it+1);
555:   VecView(x,viewer);
556:   return(0);
557: }

561: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp)
562: {
564:   PetscInt       restart;
565:   PetscReal      haptol;
566:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
567:   PetscBool      flg;

570:   PetscOptionsHead("KSP GMRES Options");
571:   PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
572:   if (flg) { KSPGMRESSetRestart(ksp,restart); }
573:   PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
574:   if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
575:   flg  = PETSC_FALSE;
576:   PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,NULL);
577:   if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
578:   PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
579:   if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
580:   PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
581:   if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
582:   PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
583:                           KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
584:   flg  = PETSC_FALSE;
585:   PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,NULL);
586:   if (flg) {
587:     PetscViewers viewers;
588:     PetscViewersCreate(PetscObjectComm((PetscObject)ksp),&viewers);
589:     KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
590:   }
591:   PetscOptionsTail();
592:   return(0);
593: }

595: extern PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal*,PetscReal*);
596: extern PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*);

600: PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
601: {
602:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

605:   if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
606:   gmres->haptol = tol;
607:   return(0);
608: }

612: PetscErrorCode  KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
613: {
614:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

617:   *max_k = gmres->max_k;
618:   return(0);
619: }

623: PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
624: {
625:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

629:   if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
630:   if (!ksp->setupstage) {
631:     gmres->max_k = max_k;
632:   } else if (gmres->max_k != max_k) {
633:     gmres->max_k    = max_k;
634:     ksp->setupstage = KSP_SETUP_NEW;
635:     /* free the data structures, then create them again */
636:     KSPReset_GMRES(ksp);
637:   }
638:   return(0);
639: }

643: PetscErrorCode  KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
644: {
646:   ((KSP_GMRES*)ksp->data)->orthog = fcn;
647:   return(0);
648: }

652: PetscErrorCode  KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
653: {
655:   *fcn = ((KSP_GMRES*)ksp->data)->orthog;
656:   return(0);
657: }

661: PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
662: {
663:   KSP_GMRES *gmres;

666:   gmres = (KSP_GMRES*)ksp->data;
667:   gmres->q_preallocate = 1;
668:   return(0);
669: }

673: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
674: {
675:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

678:   gmres->cgstype = type;
679:   return(0);
680: }

684: PetscErrorCode  KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
685: {
686:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

689:   *type = gmres->cgstype;
690:   return(0);
691: }

695: /*@
696:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
697:          in the classical Gram Schmidt orthogonalization.

699:    Logically Collective on KSP

701:    Input Parameters:
702: +  ksp - the Krylov space context
703: -  type - the type of refinement

705:   Options Database:
706: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

708:    Level: intermediate

710: .keywords: KSP, GMRES, iterative refinement

712: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
713:           KSPGMRESGetOrthogonalization()
714: @*/
715: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
716: {

722:   PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
723:   return(0);
724: }

728: /*@
729:    KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
730:          in the classical Gram Schmidt orthogonalization.

732:    Not Collective

734:    Input Parameter:
735: .  ksp - the Krylov space context

737:    Output Parameter:
738: .  type - the type of refinement

740:   Options Database:
741: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

743:    Level: intermediate

745: .keywords: KSP, GMRES, iterative refinement

747: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
748:           KSPGMRESGetOrthogonalization()
749: @*/
750: PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
751: {

756:   PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
757:   return(0);
758: }


763: /*@
764:    KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.

766:    Logically Collective on KSP

768:    Input Parameters:
769: +  ksp - the Krylov space context
770: -  restart - integer restart value

772:   Options Database:
773: .  -ksp_gmres_restart <positive integer>

775:     Note: The default value is 30.

777:    Level: intermediate

779: .keywords: KSP, GMRES, restart, iterations

781: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
782: @*/
783: PetscErrorCode  KSPGMRESSetRestart(KSP ksp, PetscInt restart)
784: {


790:   PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
791:   return(0);
792: }

796: /*@
797:    KSPGMRESGetRestart - Gets number of iterations at which GMRES, FGMRES and LGMRES restarts.

799:    Not Collective

801:    Input Parameter:
802: .  ksp - the Krylov space context

804:    Output Parameter:
805: .   restart - integer restart value

807:     Note: The default value is 30.

809:    Level: intermediate

811: .keywords: KSP, GMRES, restart, iterations

813: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
814: @*/
815: PetscErrorCode  KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
816: {

820:   PetscTryMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
821:   return(0);
822: }

826: /*@
827:    KSPGMRESSetHapTol - Sets tolerance for determining happy breakdown in GMRES, FGMRES and LGMRES.

829:    Logically Collective on KSP

831:    Input Parameters:
832: +  ksp - the Krylov space context
833: -  tol - the tolerance

835:   Options Database:
836: .  -ksp_gmres_haptol <positive real value>

838:    Note: Happy breakdown is the rare case in GMRES where an 'exact' solution is obtained after
839:          a certain number of iterations. If you attempt more iterations after this point unstable
840:          things can happen hence very occasionally you may need to set this value to detect this condition

842:    Level: intermediate

844: .keywords: KSP, GMRES, tolerance

846: .seealso: KSPSetTolerances()
847: @*/
848: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
849: {

854:   PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
855:   return(0);
856: }

858: /*MC
859:      KSPGMRES - Implements the Generalized Minimal Residual method.
860:                 (Saad and Schultz, 1986) with restart


863:    Options Database Keys:
864: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
865: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
866: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
867:                              vectors are allocated as needed)
868: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
869: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
870: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
871:                                    stability of the classical Gram-Schmidt  orthogonalization.
872: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

874:    Level: beginner

876:    Notes: Left and right preconditioning are supported, but not symmetric preconditioning.

878:    References:
879:      GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS. YOUCEF SAAD AND MARTIN H. SCHULTZ,
880:           SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986, pp. 856--869.

882: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
883:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
884:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
885:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()

887: M*/

891: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
892: {
893:   KSP_GMRES      *gmres;

897:   PetscNewLog(ksp,&gmres);
898:   ksp->data = (void*)gmres;

900:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
901:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);

903:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;
904:   ksp->ops->setup                        = KSPSetUp_GMRES;
905:   ksp->ops->solve                        = KSPSolve_GMRES;
906:   ksp->ops->reset                        = KSPReset_GMRES;
907:   ksp->ops->destroy                      = KSPDestroy_GMRES;
908:   ksp->ops->view                         = KSPView_GMRES;
909:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
910:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
911:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

913:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
914:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
915:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
916:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
917:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
918:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
919:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
920:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

922:   gmres->haptol         = 1.0e-30;
923:   gmres->q_preallocate  = 0;
924:   gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
925:   gmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
926:   gmres->nrs            = 0;
927:   gmres->sol_temp       = 0;
928:   gmres->max_k          = GMRES_DEFAULT_MAXK;
929:   gmres->Rsvd           = 0;
930:   gmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
931:   gmres->orthogwork     = 0;
932:   return(0);
933: }