Actual source code: gmres.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     This file implements GMRES (a Generalized Minimal Residual) method.  
  5:     Reference:  Saad and Schultz, 1986.


  8:     Some comments on left vs. right preconditioning, and restarts.
  9:     Left and right preconditioning.
 10:     If right preconditioning is chosen, then the problem being solved
 11:     by gmres is actually
 12:        My =  AB^-1 y = f
 13:     so the initial residual is 
 14:           r = f - Mx
 15:     Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
 16:     residual is
 17:           r = f - A x
 18:     The final solution is then
 19:           x = B^-1 y 

 21:     If left preconditioning is chosen, then the problem being solved is
 22:        My = B^-1 A x = B^-1 f,
 23:     and the initial residual is
 24:        r  = B^-1(f - Ax)

 26:     Restarts:  Restarts are basically solves with x0 not equal to zero.
 27:     Note that we can eliminate an extra application of B^-1 between
 28:     restarts as long as we don't require that the solution at the end
 29:     of an unsuccessful gmres iteration always be the solution x.
 30:  */

 32:  #include src/ksp/ksp/impls/gmres/gmresp.h
 33: #define GMRES_DELTA_DIRECTIONS 10
 34: #define GMRES_DEFAULT_MAXK     30
 35: static PetscErrorCode    GMRESGetNewVectors(KSP,PetscInt);
 36: static PetscErrorCode    GMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal*);
 37: static PetscErrorCode    BuildGmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 41: PetscErrorCode    KSPSetUp_GMRES(KSP ksp)
 42: {
 43:   PetscInt       size,hh,hes,rs,cc;
 45:   PetscInt       max_k,k;
 46:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;

 49:   if (ksp->pc_side == PC_SYMMETRIC) {
 50:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPGMRES");
 51:   }

 53:   max_k         = gmres->max_k;  /* restart size */
 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:   if (gmres->q_preallocate) {
 85:     gmres->vv_allocated   = VEC_OFFSET + 2 + max_k;
 86:     KSPGetVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,PETSC_NULL);
 87:     PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
 88:     gmres->mwork_alloc[0] = gmres->vv_allocated;
 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:   } else {
 94:     gmres->vv_allocated    = 5;
 95:     KSPGetVecs(ksp,5,&gmres->user_work[0],0,PETSC_NULL);
 96:     PetscLogObjectParents(ksp,5,gmres->user_work[0]);
 97:     gmres->mwork_alloc[0]  = 5;
 98:     gmres->nwork_alloc     = 1;
 99:     for (k=0; k<gmres->vv_allocated; k++) {
100:       gmres->vecs[k] = gmres->user_work[0][k];
101:     }
102:   }
103:   return(0);
104: }

106: /*
107:     Run gmres, possibly with restart.  Return residual history if requested.
108:     input parameters:

110: .        gmres  - structure containing parameters and work areas

112:     output parameters:
113: .        nres    - residuals (from preconditioned system) at each step.
114:                   If restarting, consider passing nres+it.  If null, 
115:                   ignored
116: .        itcount - number of iterations used.  nres[0] to nres[itcount]
117:                   are defined.  If null, ignored.
118:                   
119:     Notes:
120:     On entry, the value in vector VEC_VV(0) should be the initial residual
121:     (this allows shortcuts where the initial preconditioned residual is 0).
122:  */
125: PetscErrorCode GMREScycle(PetscInt *itcount,KSP ksp)
126: {
127:   KSP_GMRES      *gmres = (KSP_GMRES *)(ksp->data);
128:   PetscReal      res_norm,res,hapbnd,tt;
130:   PetscInt       it = 0, max_k = gmres->max_k;
131:   PetscTruth     hapend = PETSC_FALSE;

134:   VecNormalize(VEC_VV(0),&res_norm);
135:   res     = res_norm;
136:   *GRS(0) = res_norm;

138:   /* check for the convergence */
139:   PetscObjectTakeAccess(ksp);
140:   ksp->rnorm = res;
141:   PetscObjectGrantAccess(ksp);
142:   gmres->it = (it - 1);
143:   KSPLogResidualHistory(ksp,res);
144:   if (!res) {
145:     if (itcount) *itcount = 0;
146:     ksp->reason = KSP_CONVERGED_ATOL;
147:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
148:     return(0);
149:   }

151:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
152:   while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
153:     KSPLogResidualHistory(ksp,res);
154:     gmres->it = (it - 1);
155:     KSPMonitor(ksp,ksp->its,res);
156:     if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
157:       GMRESGetNewVectors(ksp,it+1);
158:     }
159:     KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);

161:     /* update hessenberg matrix and do Gram-Schmidt */
162:     (*gmres->orthog)(ksp,it);

164:     /* vv(i+1) . vv(i+1) */
165:     VecNormalize(VEC_VV(it+1),&tt);
166:     /* save the magnitude */
167:     *HH(it+1,it)    = tt;
168:     *HES(it+1,it)   = tt;

170:     /* check for the happy breakdown */
171:     hapbnd  = PetscAbsScalar(tt / *GRS(it));
172:     if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
173:     if (tt < hapbnd) {
174:       PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
175:       hapend = PETSC_TRUE;
176:     }
177:     GMRESUpdateHessenberg(ksp,it,hapend,&res);
178:     if (ksp->reason) break;

180:     it++;
181:     gmres->it  = (it-1);  /* For converged */
182:     PetscObjectTakeAccess(ksp);
183:     ksp->its++;
184:     ksp->rnorm = res;
185:     PetscObjectGrantAccess(ksp);

187:     (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);

189:     /* Catch error in happy breakdown and signal convergence and break from loop */
190:     if (hapend) {
191:       if (!ksp->reason) {
192:         SETERRQ1(0,"You reached the happy break down, but convergence was not indicated. Residual norm = %G",res);
193:       }
194:       break;
195:     }
196:   }

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

204:   if (itcount) *itcount    = it;


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

215:   return(0);
216: }

220: PetscErrorCode KSPSolve_GMRES(KSP ksp)
221: {
223:   PetscInt       its,itcount;
224:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;
225:   PetscTruth     guess_zero = ksp->guess_zero;

228:   if (ksp->calc_sings && !gmres->Rsvd) {
229:     SETERRQ(PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
230:   }

232:   PetscObjectTakeAccess(ksp);
233:   ksp->its = 0;
234:   PetscObjectGrantAccess(ksp);

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

254: PetscErrorCode KSPDestroy_GMRES_Internal(KSP ksp)
255: {
256:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
258:   PetscInt       i;

261:   /* Free the Hessenberg matrix */
262:   PetscFree(gmres->hh_origin);

264:   /* Free the pointer to user variables */
265:   PetscFree(gmres->vecs);

267:   /* free work vectors */
268:   for (i=0; i<gmres->nwork_alloc; i++) {
269:     VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
270:   }
271:   PetscFree(gmres->user_work);
272:   PetscFree(gmres->mwork_alloc);
273:   PetscFree(gmres->nrs);
274:   if (gmres->sol_temp) {
275:     VecDestroy(gmres->sol_temp);
276:   }
277:   PetscFree(gmres->Rsvd);
278:   PetscFree(gmres->Dsvd);
279:   PetscFree(gmres->orthogwork);
280:   gmres->sol_temp       = 0;
281:   gmres->vv_allocated   = 0;
282:   gmres->vecs_allocated = 0;
283:   gmres->sol_temp       = 0;
284:   return(0);
285: }

289: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
290: {
291:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

295:   KSPDestroy_GMRES_Internal(ksp);
296:   PetscFree(gmres);
297:   return(0);
298: }
299: /*
300:     BuildGmresSoln - create the solution from the starting vector and the
301:     current iterates.

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

309:      This is an internal routine that knows about the GMRES internals.
310:  */
313: static PetscErrorCode BuildGmresSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
314: {
315:   PetscScalar    tt;
317:   PetscInt       ii,k,j;
318:   KSP_GMRES      *gmres = (KSP_GMRES *)(ksp->data);

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

323:   /* If it is < 0, no gmres steps have been performed */
324:   if (it < 0) {
325:     if (vdest != vs) {
326:       VecCopy(vs,vdest);
327:     }
328:     return(0);
329:   }
330:   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)));
331:   if (*HH(it,it) != 0.0) {
332:     nrs[it] = *GRS(it) / *HH(it,it);
333:   } else {
334:     nrs[it] = 0.0;
335:   }
336:   for (ii=1; ii<=it; ii++) {
337:     k   = it - ii;
338:     tt  = *GRS(k);
339:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
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.
357:  */
360: static PetscErrorCode GMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth 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: static PetscErrorCode GMRESGetNewVectors(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:   BuildGmresSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
476:   *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:   PetscTruth     iascii,isstring;

490:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
491:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&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(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(PETSC_ERR_SUP,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
518:   }
519:   return(0);
520: }

524: /*@C
525:    KSPGMRESKrylovMonitor - 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: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
541: @*/
542: PetscErrorCode  KSPGMRESKrylovMonitor(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;

551:   PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
552:   PetscViewerSetType(viewer,PETSC_VIEWER_DRAW);

554:   x      = VEC_VV(gmres->it+1);
555:   VecView(x,viewer);

557:   return(0);
558: }

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

571:   PetscOptionsHead("KSP GMRES Options");
572:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
573:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
574:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
575:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
576:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
577:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
578:     PetscOptionsTruthGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
579:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
580:     PetscOptionsTruthGroupEnd("-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:     PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
585:     if (flg) {
586:       PetscViewers viewers;
587:       PetscViewersCreate(ksp->comm,&viewers);
588:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(PetscErrorCode (*)(void*))PetscViewersDestroy);
589:     }
590:   PetscOptionsTail();
591:   return(0);
592: }

594: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
595: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);


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

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

615: PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
616: {
617:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;

621:   if (max_k < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
622:   if (!ksp->setupcalled) {
623:     gmres->max_k = max_k;
624:   } else if (gmres->max_k != max_k) {
625:      gmres->max_k = max_k;
626:      ksp->setupcalled = 0;
627:      /* free the data structures, then create them again */
628:      KSPDestroy_GMRES_Internal(ksp);
629:   }
630:   return(0);
631: }

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

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: }

664: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
665: {
666:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

669:   gmres->cgstype = type;
670:   return(0);
671: }

676: /*@
677:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
678:          in the classical Gram Schmidt orthogonalization.
679:    of the preconditioned problem.

681:    Collective on KSP

683:    Input Parameters:
684: +  ksp - the Krylov space context
685: -  type - the type of refinement

687:   Options Database:
688: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

690:    Level: intermediate

692: .keywords: KSP, GMRES, iterative refinement

694: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization()
695: @*/
696: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
697: {
698:   PetscErrorCode ierr,(*f)(KSP,KSPGMRESCGSRefinementType);

702:   PetscObjectQueryFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",(void (**)(void))&f);
703:   if (f) {
704:     (*f)(ksp,type);
705:   }
706:   return(0);
707: }

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

714:    Collective on KSP

716:    Input Parameters:
717: +  ksp - the Krylov space context
718: -  restart - integer restart value

720:   Options Database:
721: .  -ksp_gmres_restart <positive integer>

723:     Note: The default value is 30.

725:    Level: intermediate

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

729: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors()
730: @*/
731: PetscErrorCode  KSPGMRESSetRestart(KSP ksp, PetscInt restart)
732: {

736:   PetscTryMethod(ksp,KSPGMRESSetRestart_C,(KSP,PetscInt),(ksp,restart));
737:   return(0);
738: }

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

745:    Collective on KSP

747:    Input Parameters:
748: +  ksp - the Krylov space context
749: -  tol - the tolerance

751:   Options Database:
752: .  -ksp_gmres_haptol <positive real value>

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

758:    Level: intermediate

760: .keywords: KSP, GMRES, tolerance

762: .seealso: KSPSetTolerances()
763: @*/
764: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
765: {

769:   PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)));
770:   return(0);
771: }

773: /*MC
774:      KSPGMRES - Implements the Generalized Minimal Residual method.  
775:                 (Saad and Schultz, 1986) with restart


778:    Options Database Keys:
779: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
780: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
781: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 
782:                              vectors are allocated as needed)
783: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
784: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
785: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 
786:                                    stability of the classical Gram-Schmidt  orthogonalization.
787: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

789:    Level: beginner


792: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
793:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
794:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
795:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor()

797: M*/

802: PetscErrorCode  KSPCreate_GMRES(KSP ksp)
803: {
804:   KSP_GMRES      *gmres;

808:   PetscNew(KSP_GMRES,&gmres);
809:   PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
810:   ksp->data                              = (void*)gmres;
811:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;

813:   ksp->ops->setup                        = KSPSetUp_GMRES;
814:   ksp->ops->solve                        = KSPSolve_GMRES;
815:   ksp->ops->destroy                      = KSPDestroy_GMRES;
816:   ksp->ops->view                         = KSPView_GMRES;
817:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
818:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
819:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

821:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
822:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
823:                                      KSPGMRESSetPreAllocateVectors_GMRES);
824:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
825:                                     "KSPGMRESSetOrthogonalization_GMRES",
826:                                      KSPGMRESSetOrthogonalization_GMRES);
827:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
828:                                     "KSPGMRESSetRestart_GMRES",
829:                                      KSPGMRESSetRestart_GMRES);
830:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
831:                                     "KSPGMRESSetHapTol_GMRES",
832:                                      KSPGMRESSetHapTol_GMRES);
833:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
834:                                     "KSPGMRESSetCGSRefinementType_GMRES",
835:                                      KSPGMRESSetCGSRefinementType_GMRES);

837:   gmres->haptol              = 1.0e-30;
838:   gmres->q_preallocate       = 0;
839:   gmres->delta_allocate      = GMRES_DELTA_DIRECTIONS;
840:   gmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
841:   gmres->nrs                 = 0;
842:   gmres->sol_temp            = 0;
843:   gmres->max_k               = GMRES_DEFAULT_MAXK;
844:   gmres->Rsvd                = 0;
845:   gmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
846:   gmres->orthogwork          = 0;
847:   return(0);
848: }