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:   }
231:   if (ksp->normtype != KSP_NORM_PRECONDITIONED) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Currently can use GMRES with only preconditioned residual (right preconditioning not coded)");

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

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

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

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

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

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

290: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
291: {

295:   KSPDestroy_GMRES_Internal(ksp);
296:   PetscFree(ksp->data);
297:   /* clear composed functions */
298:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C","",PETSC_NULL);
299:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C","",PETSC_NULL);
300:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C","",PETSC_NULL);
301:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C","",PETSC_NULL);
302:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C","",PETSC_NULL);
303:   return(0);
304: }
305: /*
306:     BuildGmresSoln - create the solution from the starting vector and the
307:     current iterates.

309:     Input parameters:
310:         nrs - work area of size it + 1.
311:         vs  - index of initial guess
312:         vdest - index of result.  Note that vs may == vdest (replace
313:                 guess with the solution).

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

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

329:   /* If it is < 0, no gmres steps have been performed */
330:   if (it < 0) {
331:     if (vdest != vs) {
332:       VecCopy(vs,vdest);
333:     }
334:     return(0);
335:   }
336:   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)));
337:   if (*HH(it,it) != 0.0) {
338:     nrs[it] = *GRS(it) / *HH(it,it);
339:   } else {
340:     nrs[it] = 0.0;
341:   }
342:   for (ii=1; ii<=it; ii++) {
343:     k   = it - ii;
344:     tt  = *GRS(k);
345:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
346:     nrs[k]   = tt / *HH(k,k);
347:   }

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

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

373:   hh  = HH(0,it);
374:   cc  = CC(0);
375:   ss  = SS(0);

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

390:   /*
391:     compute the new plane rotation, and apply it to:
392:      1) the right-hand-side of the Hessenberg system
393:      2) the new column of the Hessenberg matrix
394:     thus obtaining the updated value of the residual
395:   */
396:   if (!hapend) {
397: #if defined(PETSC_USE_COMPLEX)
398:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
399: #else
400:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
401: #endif
402:     if (tt == 0.0) {
403:       ksp->reason = KSP_DIVERGED_NULL;
404:       return(0);
405:     }
406:     *cc       = *hh / tt;
407:     *ss       = *(hh+1) / tt;
408:     *GRS(it+1) = - (*ss * *GRS(it));
409: #if defined(PETSC_USE_COMPLEX)
410:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
411:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);
412: #else
413:     *GRS(it)   = *cc * *GRS(it);
414:     *hh       = *cc * *hh + *ss * *(hh+1);
415: #endif
416:     *res      = PetscAbsScalar(*GRS(it+1));
417:   } else {
418:     /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
419:             another rotation matrix (so RH doesn't change).  The new residual is 
420:             always the new sine term times the residual from last time (GRS(it)), 
421:             but now the new sine rotation would be zero...so the residual should
422:             be zero...so we will multiply "zero" by the last residual.  This might
423:             not be exactly what we want to do here -could just return "zero". */
424: 
425:     *res = 0.0;
426:   }
427:   return(0);
428: }
429: /*
430:    This routine allocates more work vectors, starting from VEC_VV(it).
431:  */
434: static PetscErrorCode GMRESGetNewVectors(KSP ksp,PetscInt it)
435: {
436:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;
438:   PetscInt       nwork = gmres->nwork_alloc,k,nalloc;

441:   nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
442:   /* Adjust the number to allocate to make sure that we don't exceed the
443:     number of available slots */
444:   if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
445:     nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
446:   }
447:   if (!nalloc) return(0);

449:   gmres->vv_allocated += nalloc;
450:   KSPGetVecs(ksp,nalloc,&gmres->user_work[nwork],0,PETSC_NULL);
451:   PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
452:   gmres->mwork_alloc[nwork] = nalloc;
453:   for (k=0; k<nalloc; k++) {
454:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
455:   }
456:   gmres->nwork_alloc++;
457:   return(0);
458: }

462: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec  ptr,Vec *result)
463: {
464:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;

468:   if (!ptr) {
469:     if (!gmres->sol_temp) {
470:       VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
471:       PetscLogObjectParent(ksp,gmres->sol_temp);
472:     }
473:     ptr = gmres->sol_temp;
474:   }
475:   if (!gmres->nrs) {
476:     /* allocate the work area */
477:     PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
478:     PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
479:   }

481:   BuildGmresSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
482:   *result = ptr;
483:   return(0);
484: }

488: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
489: {
490:   KSP_GMRES      *gmres = (KSP_GMRES *)ksp->data;
491:   const char     *cstr;
493:   PetscTruth     iascii,isstring;

496:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
497:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
498:   if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
499:     switch (gmres->cgstype) {
500:       case (KSP_GMRES_CGS_REFINE_NEVER):
501:         cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
502:         break;
503:       case (KSP_GMRES_CGS_REFINE_ALWAYS):
504:         cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
505:         break;
506:       case (KSP_GMRES_CGS_REFINE_IFNEEDED):
507:         cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
508:         break;
509:       default:
510:         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
511:     }
512:   } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
513:     cstr = "Modified Gram-Schmidt Orthogonalization";
514:   } else {
515:     cstr = "unknown orthogonalization";
516:   }
517:   if (iascii) {
518:     PetscViewerASCIIPrintf(viewer,"  GMRES: restart=%D, using %s\n",gmres->max_k,cstr);
519:     PetscViewerASCIIPrintf(viewer,"  GMRES: happy breakdown tolerance %G\n",gmres->haptol);
520:   } else if (isstring) {
521:     PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
522:   } else {
523:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
524:   }
525:   return(0);
526: }

530: /*@C
531:    KSPGMRESMonitorKrylov - Calls VecView() for each direction in the 
532:    GMRES accumulated Krylov space.

534:    Collective on KSP

536:    Input Parameters:
537: +  ksp - the KSP context
538: .  its - iteration number
539: .  fgnorm - 2-norm of residual (or gradient)
540: -  a viewers object created with PetscViewersCreate()

542:    Level: intermediate

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

546: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
547: @*/
548: PetscErrorCode  KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
549: {
550:   PetscViewers   viewers = (PetscViewers)dummy;
551:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
553:   Vec            x;
554:   PetscViewer    viewer;

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

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

563:   return(0);
564: }

568: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp)
569: {
571:   PetscInt       restart;
572:   PetscReal      haptol;
573:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
574:   PetscTruth     flg;

577:   PetscOptionsHead("KSP GMRES Options");
578:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
579:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
580:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
581:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
582:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
583:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
584:     PetscOptionsTruthGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
585:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
586:     PetscOptionsTruthGroupEnd("-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:     PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",&flg);
591:     if (flg) {
592:       PetscViewers viewers;
593:       PetscViewersCreate(ksp->comm,&viewers);
594:       KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void*))PetscViewersDestroy);
595:     }
596:   PetscOptionsTail();
597:   return(0);
598: }

600: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
601: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);


607: PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
608: {
609:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;

612:   if (tol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
613:   gmres->haptol = tol;
614:   return(0);
615: }

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

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

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

656: PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
657: {
658:   KSP_GMRES *gmres;

661:   gmres = (KSP_GMRES *)ksp->data;
662:   gmres->q_preallocate = 1;
663:   return(0);
664: }

670: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
671: {
672:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

675:   gmres->cgstype = type;
676:   return(0);
677: }

682: /*@
683:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
684:          in the classical Gram Schmidt orthogonalization.
685:    of the preconditioned problem.

687:    Collective on KSP

689:    Input Parameters:
690: +  ksp - the Krylov space context
691: -  type - the type of refinement

693:   Options Database:
694: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

696:    Level: intermediate

698: .keywords: KSP, GMRES, iterative refinement

700: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization()
701: @*/
702: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
703: {
704:   PetscErrorCode ierr,(*f)(KSP,KSPGMRESCGSRefinementType);

708:   PetscObjectQueryFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",(void (**)(void))&f);
709:   if (f) {
710:     (*f)(ksp,type);
711:   }
712:   return(0);
713: }

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

720:    Collective on KSP

722:    Input Parameters:
723: +  ksp - the Krylov space context
724: -  restart - integer restart value

726:   Options Database:
727: .  -ksp_gmres_restart <positive integer>

729:     Note: The default value is 30.

731:    Level: intermediate

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

735: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors()
736: @*/
737: PetscErrorCode  KSPGMRESSetRestart(KSP ksp, PetscInt restart)
738: {

742:   PetscTryMethod(ksp,KSPGMRESSetRestart_C,(KSP,PetscInt),(ksp,restart));
743:   return(0);
744: }

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

751:    Collective on KSP

753:    Input Parameters:
754: +  ksp - the Krylov space context
755: -  tol - the tolerance

757:   Options Database:
758: .  -ksp_gmres_haptol <positive real value>

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

764:    Level: intermediate

766: .keywords: KSP, GMRES, tolerance

768: .seealso: KSPSetTolerances()
769: @*/
770: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
771: {

775:   PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)));
776:   return(0);
777: }

779: /*MC
780:      KSPGMRES - Implements the Generalized Minimal Residual method.  
781:                 (Saad and Schultz, 1986) with restart


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

795:    Level: beginner


798: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
799:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
800:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
801:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov()

803: M*/

808: PetscErrorCode  KSPCreate_GMRES(KSP ksp)
809: {
810:   KSP_GMRES      *gmres;

814:   PetscNew(KSP_GMRES,&gmres);
815:   PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
816:   ksp->data                              = (void*)gmres;


819:   ksp->normtype                          = KSP_NORM_PRECONDITIONED;
820:   ksp->pc_side                           = PC_LEFT;

822:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;
823:   ksp->ops->setup                        = KSPSetUp_GMRES;
824:   ksp->ops->solve                        = KSPSolve_GMRES;
825:   ksp->ops->destroy                      = KSPDestroy_GMRES;
826:   ksp->ops->view                         = KSPView_GMRES;
827:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
828:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
829:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

831:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
832:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
833:                                      KSPGMRESSetPreAllocateVectors_GMRES);
834:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
835:                                     "KSPGMRESSetOrthogonalization_GMRES",
836:                                      KSPGMRESSetOrthogonalization_GMRES);
837:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
838:                                     "KSPGMRESSetRestart_GMRES",
839:                                      KSPGMRESSetRestart_GMRES);
840:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
841:                                     "KSPGMRESSetHapTol_GMRES",
842:                                      KSPGMRESSetHapTol_GMRES);
843:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
844:                                     "KSPGMRESSetCGSRefinementType_GMRES",
845:                                      KSPGMRESSetCGSRefinementType_GMRES);

847:   gmres->haptol              = 1.0e-30;
848:   gmres->q_preallocate       = 0;
849:   gmres->delta_allocate      = GMRES_DELTA_DIRECTIONS;
850:   gmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
851:   gmres->nrs                 = 0;
852:   gmres->sol_temp            = 0;
853:   gmres->max_k               = GMRES_DEFAULT_MAXK;
854:   gmres->Rsvd                = 0;
855:   gmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
856:   gmres->orthogwork          = 0;
857:   return(0);
858: }