Actual source code: gmres.c

  1: /*$Id: gmres.c,v 1.176 2001/08/07 03:03:51 balay Exp $*/

  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 int    GMRESGetNewVectors(KSP,int);
 36: static int    GMRESUpdateHessenberg(KSP,int,PetscTruth,PetscReal*);
 37: static int    BuildGmresSoln(PetscScalar*,Vec,Vec,KSP,int);

 41: int    KSPSetUp_GMRES(KSP ksp)
 42: {
 43:   unsigned  int size,hh,hes,rs,cc;
 44:   int       ierr,max_k,k;
 45:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;

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

 52:   max_k         = gmres->max_k;
 53:   hh            = (max_k + 2) * (max_k + 1);
 54:   hes           = (max_k + 1) * (max_k + 1);
 55:   rs            = (max_k + 2);
 56:   cc            = (max_k + 1);
 57:   size          = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);

 59:   PetscMalloc(size,&gmres->hh_origin);
 60:   PetscMemzero(gmres->hh_origin,size);
 61:   PetscLogObjectMemory(ksp,size);
 62:   gmres->hes_origin = gmres->hh_origin + hh;
 63:   gmres->rs_origin  = gmres->hes_origin + hes;
 64:   gmres->cc_origin  = gmres->rs_origin + rs;
 65:   gmres->ss_origin  = gmres->cc_origin + cc;

 67:   if (ksp->calc_sings) {
 68:     /* Allocate workspace to hold Hessenberg matrix needed by lapack */
 69:     size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
 70:     PetscMalloc(size,&gmres->Rsvd);
 71:     PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
 72:     PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
 73:   }

 75:   /* Allocate array to hold pointers to user vectors.  Note that we need
 76:    4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
 77:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&gmres->vecs);
 78:   gmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
 79:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&gmres->user_work);
 80:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(int),&gmres->mwork_alloc);
 81:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)+sizeof(int)));

 83:   if (gmres->q_preallocate) {
 84:     gmres->vv_allocated   = VEC_OFFSET + 2 + max_k;
 85:     VecDuplicateVecs(VEC_RHS,gmres->vv_allocated,&gmres->user_work[0]);
 86:     PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
 87:     gmres->mwork_alloc[0] = gmres->vv_allocated;
 88:     gmres->nwork_alloc    = 1;
 89:     for (k=0; k<gmres->vv_allocated; k++) {
 90:       gmres->vecs[k] = gmres->user_work[0][k];
 91:     }
 92:   } else {
 93:     gmres->vv_allocated    = 5;
 94:     VecDuplicateVecs(ksp->vec_rhs,5,&gmres->user_work[0]);
 95:     PetscLogObjectParents(ksp,5,gmres->user_work[0]);
 96:     gmres->mwork_alloc[0]  = 5;
 97:     gmres->nwork_alloc     = 1;
 98:     for (k=0; k<gmres->vv_allocated; k++) {
 99:       gmres->vecs[k] = gmres->user_work[0][k];
100:     }
101:   }
102:   return(0);
103: }

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

109: .        gmres  - structure containing parameters and work areas

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

132:   VecNormalize(VEC_VV(0),&res_norm);
133:   res     = res_norm;
134:   *GRS(0) = res_norm;

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

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

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

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

168:     /* check for the happy breakdown */
169:     hapbnd  = PetscAbsScalar(tt / *GRS(it));
170:     if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
171:     if (tt < hapbnd) {
172:       PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",hapbnd,tt);
173:       hapend = PETSC_TRUE;
174:     }
175:     GMRESUpdateHessenberg(ksp,it,hapend,&res);
176:     it++;
177:     gmres->it  = (it-1);  /* For converged */
178:     PetscObjectTakeAccess(ksp);
179:     ksp->its++;
180:     ksp->rnorm = res;
181:     PetscObjectGrantAccess(ksp);

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

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

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

200:   if (itcount) *itcount    = it;


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

211:   return(0);
212: }

216: int KSPSolve_GMRES(KSP ksp)
217: {
218:   int        ierr,its,itcount;
219:   KSP_GMRES  *gmres = (KSP_GMRES *)ksp->data;
220:   PetscTruth guess_zero = ksp->guess_zero;

223:   if (ksp->calc_sings && !gmres->Rsvd) {
224:     SETERRQ(1,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
225:   }

227:   PetscObjectTakeAccess(ksp);
228:   ksp->its = 0;
229:   PetscObjectGrantAccess(ksp);

231:   itcount     = 0;
232:   ksp->reason = KSP_CONVERGED_ITERATING;
233:   while (!ksp->reason) {
234:     KSPInitialResidual(ksp,VEC_SOLN,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),VEC_RHS);
235:     GMREScycle(&its,ksp);
236:     itcount += its;
237:     if (itcount >= ksp->max_it) {
238:       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: int KSPDestroy_GMRES_Internal(KSP ksp)
250: {
251:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
252:   int       i,ierr;

255:   /* Free the Hessenberg matrix */
256:   if (gmres->hh_origin) {PetscFree(gmres->hh_origin);}

258:   /* Free the pointer to user variables */
259:   if (gmres->vecs) {PetscFree(gmres->vecs);}

261:   /* free work vectors */
262:   for (i=0; i<gmres->nwork_alloc; i++) {
263:     VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
264:   }
265:   if (gmres->user_work)  {PetscFree(gmres->user_work);}
266:   if (gmres->mwork_alloc) {PetscFree(gmres->mwork_alloc);}
267:   if (gmres->nrs) {PetscFree(gmres->nrs);}
268:   if (gmres->sol_temp) {VecDestroy(gmres->sol_temp);}
269:   if (gmres->Rsvd) {PetscFree(gmres->Rsvd);}
270:   if (gmres->Dsvd) {PetscFree(gmres->Dsvd);}

272:   return(0);
273: }

277: int KSPDestroy_GMRES(KSP ksp)
278: {
279:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
280:   int       ierr;

283:   KSPDestroy_GMRES_Internal(ksp);
284:   PetscFree(gmres);
285:   return(0);
286: }
287: /*
288:     BuildGmresSoln - create the solution from the starting vector and the
289:     current iterates.

291:     Input parameters:
292:         nrs - work area of size it + 1.
293:         vs  - index of initial guess
294:         vdest - index of result.  Note that vs may == vdest (replace
295:                 guess with the solution).

297:      This is an internal routine that knows about the GMRES internals.
298:  */
301: static int BuildGmresSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,int it)
302: {
303:   PetscScalar tt,zero = 0.0,one = 1.0;
304:   int         ierr,ii,k,j;
305:   KSP_GMRES   *gmres = (KSP_GMRES *)(ksp->data);

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

310:   /* If it is < 0, no gmres steps have been performed */
311:   if (it < 0) {
312:     if (vdest != vs) {
313:       VecCopy(vs,vdest);
314:     }
315:     return(0);
316:   }
317:   if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
318:   if (*HH(it,it) != 0.0) {
319:     nrs[it] = *GRS(it) / *HH(it,it);
320:   } else {
321:     nrs[it] = 0.0;
322:   }
323:   for (ii=1; ii<=it; ii++) {
324:     k   = it - ii;
325:     tt  = *GRS(k);
326:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
327:     nrs[k]   = tt / *HH(k,k);
328:   }

330:   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
331:   VecSet(&zero,VEC_TEMP);
332:   VecMAXPY(it+1,nrs,VEC_TEMP,&VEC_VV(0));

334:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
335:   /* add solution to previous solution */
336:   if (vdest != vs) {
337:     VecCopy(vs,vdest);
338:   }
339:   VecAXPY(&one,VEC_TEMP,vdest);
340:   return(0);
341: }
342: /*
343:    Do the scalar work for the orthogonalization.  Return new residual.
344:  */
347: static int GMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
348: {
349:   PetscScalar *hh,*cc,*ss,tt;
350:   int         j;
351:   KSP_GMRES   *gmres = (KSP_GMRES *)(ksp->data);

354:   hh  = HH(0,it);
355:   cc  = CC(0);
356:   ss  = SS(0);

358:   /* Apply all the previously computed plane rotations to the new column
359:      of the Hessenberg matrix */
360:   for (j=1; j<=it; j++) {
361:     tt  = *hh;
362: #if defined(PETSC_USE_COMPLEX)
363:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
364: #else
365:     *hh = *cc * tt + *ss * *(hh+1);
366: #endif
367:     hh++;
368:     *hh = *cc++ * *hh - (*ss++ * tt);
369:   }

371:   /*
372:     compute the new plane rotation, and apply it to:
373:      1) the right-hand-side of the Hessenberg system
374:      2) the new column of the Hessenberg matrix
375:     thus obtaining the updated value of the residual
376:   */
377:   if (!hapend) {
378: #if defined(PETSC_USE_COMPLEX)
379:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
380: #else
381:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
382: #endif
383:     if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
384:     *cc       = *hh / tt;
385:     *ss       = *(hh+1) / tt;
386:     *GRS(it+1) = - (*ss * *GRS(it));
387: #if defined(PETSC_USE_COMPLEX)
388:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
389:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);
390: #else
391:     *GRS(it)   = *cc * *GRS(it);
392:     *hh       = *cc * *hh + *ss * *(hh+1);
393: #endif
394:     *res      = PetscAbsScalar(*GRS(it+1));
395:   } else {
396:     /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
397:             another rotation matrix (so RH doesn't change).  The new residual is 
398:             always the new sine term times the residual from last time (GRS(it)), 
399:             but now the new sine rotation would be zero...so the residual should
400:             be zero...so we will multiply "zero" by the last residual.  This might
401:             not be exactly what we want to do here -could just return "zero". */
402: 
403:     *res = 0.0;
404:   }
405:   return(0);
406: }
407: /*
408:    This routine allocates more work vectors, starting from VEC_VV(it).
409:  */
412: static int GMRESGetNewVectors(KSP ksp,int it)
413: {
414:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
415:   int       nwork = gmres->nwork_alloc,k,nalloc,ierr;

418:   nalloc = gmres->delta_allocate;
419:   /* Adjust the number to allocate to make sure that we don't exceed the
420:     number of available slots */
421:   if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
422:     nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
423:   }
424:   if (!nalloc) return(0);

426:   gmres->vv_allocated += nalloc;
427:   VecDuplicateVecs(ksp->vec_rhs,nalloc,&gmres->user_work[nwork]);
428:   PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
429:   gmres->mwork_alloc[nwork] = nalloc;
430:   for (k=0; k<nalloc; k++) {
431:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
432:   }
433:   gmres->nwork_alloc++;
434:   return(0);
435: }

439: int KSPBuildSolution_GMRES(KSP ksp,Vec  ptr,Vec *result)
440: {
441:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
442:   int       ierr;

445:   if (!ptr) {
446:     if (!gmres->sol_temp) {
447:       VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
448:       PetscLogObjectParent(ksp,gmres->sol_temp);
449:     }
450:     ptr = gmres->sol_temp;
451:   }
452:   if (!gmres->nrs) {
453:     /* allocate the work area */
454:     PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
455:     PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
456:   }

458:   BuildGmresSoln(gmres->nrs,VEC_SOLN,ptr,ksp,gmres->it);
459:   *result = ptr;
460:   return(0);
461: }

465: int KSPView_GMRES(KSP ksp,PetscViewer viewer)
466: {
467:   KSP_GMRES  *gmres = (KSP_GMRES *)ksp->data;
468:   const char *cstr;
469:   int        ierr;
470:   PetscTruth isascii,isstring;

473:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
474:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
475:   if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
476:     if (gmres->cgstype == KSP_GMRES_CGS_REFINE_NEVER) {
477:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
478:     } else if (gmres->cgstype == KSP_GMRES_CGS_REFINE_ALWAYS) {
479:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
480:     } else {
481:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
482:     }
483:   } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
484:     cstr = "Modified Gram-Schmidt Orthogonalization";
485:   } else {
486:     cstr = "unknown orthogonalization";
487:   }
488:   if (isascii) {
489:     PetscViewerASCIIPrintf(viewer,"  GMRES: restart=%d, using %s\n",gmres->max_k,cstr);
490:     PetscViewerASCIIPrintf(viewer,"  GMRES: happy breakdown tolerance %g\n",gmres->haptol);
491:   } else if (isstring) {
492:     PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,gmres->max_k);
493:   } else {
494:     SETERRQ1(1,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
495:   }
496:   return(0);
497: }

501: /*@C
502:    KSPGMRESKrylovMonitor - Calls VecView() for each direction in the 
503:    GMRES accumulated Krylov space.

505:    Collective on KSP

507:    Input Parameters:
508: +  ksp - the KSP context
509: .  its - iteration number
510: .  fgnorm - 2-norm of residual (or gradient)
511: -  a viewers object created with PetscViewersCreate()

513:    Level: intermediate

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

517: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
518: @*/
519: int KSPGMRESKrylovMonitor(KSP ksp,int its,PetscReal fgnorm,void *dummy)
520: {
521:   PetscViewers viewers = (PetscViewers)dummy;
522:   KSP_GMRES    *gmres = (KSP_GMRES*)ksp->data;
523:   int          ierr;
524:   Vec          x;
525:   PetscViewer  viewer;

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

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

534:   return(0);
535: }

539: int KSPSetFromOptions_GMRES(KSP ksp)
540: {
541:   int             ierr,restart,indx;
542:   PetscReal       haptol;
543:   KSP_GMRES       *gmres = (KSP_GMRES*)ksp->data;
544:   PetscTruth      flg;
545:   const char      *types[] = {"never","ifneeded","always"};

548:   PetscOptionsHead("KSP GMRES Options");
549:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
550:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
551:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
552:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
553:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
554:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
555:     PetscOptionsLogicalGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
556:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
557:     PetscOptionsLogicalGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
558:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
559:     PetscOptionsEList("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType()",types,3,types[(int)gmres->cgstype],&indx,&flg);
560:     if (flg) {
561:       KSPGMRESSetCGSRefinementType(ksp,(KSPGMRESCGSRefinementType)indx);
562:     }

564:     PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
565:     if (flg) {
566:       PetscViewers viewers;
567:       PetscViewersCreate(ksp->comm,&viewers);
568:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
569:     }
570:   PetscOptionsTail();
571:   return(0);
572: }

574: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
575: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);


578: EXTERN_C_BEGIN
581: int KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
582: {
583:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;

586:   if (tol < 0.0) SETERRQ(1,"Tolerance must be non-negative");
587:   gmres->haptol = tol;
588:   return(0);
589: }
590: EXTERN_C_END

592: EXTERN_C_BEGIN
595: int KSPGMRESSetRestart_GMRES(KSP ksp,int max_k)
596: {
597:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
598:   int       ierr;

601:   if (max_k < 1) SETERRQ(1,"Restart must be positive");
602:   if (!ksp->setupcalled) {
603:     gmres->max_k = max_k;
604:   } else if (gmres->max_k != max_k) {
605:      gmres->max_k = max_k;
606:      ksp->setupcalled = 0;
607:      /* free the data structures, then create them again */
608:      KSPDestroy_GMRES_Internal(ksp);
609:   }

611:   return(0);
612: }
613: EXTERN_C_END

615: typedef int (*FCN)(KSP,int); /* force argument to next function to not be extern C*/
616: EXTERN_C_BEGIN
619: int KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
620: {
623:   ((KSP_GMRES *)ksp->data)->orthog = fcn;
624:   return(0);
625: }
626: EXTERN_C_END

628: EXTERN_C_BEGIN
631: int KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
632: {
633:   KSP_GMRES *gmres;

636:   gmres = (KSP_GMRES *)ksp->data;
637:   gmres->q_preallocate = 1;
638:   return(0);
639: }
640: EXTERN_C_END

642: EXTERN_C_BEGIN
645: int KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
646: {
647:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

650:   gmres->cgstype = type;
651:   return(0);
652: }
653: EXTERN_C_END

657: /*@
658:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
659:          in the classical Gram Schmidt orthogonalization.
660:    of the preconditioned problem.

662:    Collective on KSP

664:    Input Parameters:
665: +  ksp - the Krylov space context
666: -  type - the type of refinement

668:   Options Database:
669: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

671:    Level: intermediate

673: .keywords: KSP, GMRES, iterative refinement

675: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization()
676: @*/
677: int KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
678: {
679:   int ierr,(*f)(KSP,KSPGMRESCGSRefinementType);

683:   PetscObjectQueryFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",(void (**)(void))&f);
684:   if (f) {
685:     (*f)(ksp,type);
686:   }
687:   return(0);
688: }

692: /*@C
693:    KSPGMRESSetRestart - Sets number of iterations at which GMRES, FGMRES and LGMRES restarts.

695:    Collective on KSP

697:    Input Parameters:
698: +  ksp - the Krylov space context
699: -  restart - integer restart value

701:   Options Database:
702: .  -ksp_gmres_restart <positive integer>

704:     Note: The default value is 30.

706:    Level: intermediate

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

710: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors()
711: @*/
712: int KSPGMRESSetRestart(KSP ksp, int restart)
713: {

717:   PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,int),((ksp),(restart)));
718:   return(0);
719: }

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

726:    Collective on KSP

728:    Input Parameters:
729: +  ksp - the Krylov space context
730: -  tol - the tolerance

732:   Options Database:
733: .  -ksp_gmres_haptol <positive real value>

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

739:    Level: intermediate

741: .keywords: KSP, GMRES, tolerance

743: .seealso: KSPSetTolerances()
744: @*/
745: int KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
746: {

750:   PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)));
751:   return(0);
752: }

754: /*MC
755:      KSPGMRES - Implements the Generalized Minimal Residual method.  
756:                 (Saad and Schultz, 1986) with restart


759:    Options Database Keys:
760: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
761: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
762: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 
763:                              vectors are allocated as needed)
764: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
765: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
766: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 
767:                                    stability of the classical Gram-Schmidt  orthogonalization.
768: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

770:    Level: beginner


773: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
774:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
775:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
776:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor()

778: M*/

780: EXTERN_C_BEGIN
783: int KSPCreate_GMRES(KSP ksp)
784: {
785:   KSP_GMRES *gmres;
786:   int       ierr;

789:   PetscNew(KSP_GMRES,&gmres);
790:   PetscMemzero(gmres,sizeof(KSP_GMRES));
791:   PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
792:   ksp->data                              = (void*)gmres;
793:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;

795:   ksp->ops->setup                        = KSPSetUp_GMRES;
796:   ksp->ops->solve                        = KSPSolve_GMRES;
797:   ksp->ops->destroy                      = KSPDestroy_GMRES;
798:   ksp->ops->view                         = KSPView_GMRES;
799:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
800:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
801:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

803:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
804:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
805:                                      KSPGMRESSetPreAllocateVectors_GMRES);
806:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
807:                                     "KSPGMRESSetOrthogonalization_GMRES",
808:                                      KSPGMRESSetOrthogonalization_GMRES);
809:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
810:                                     "KSPGMRESSetRestart_GMRES",
811:                                      KSPGMRESSetRestart_GMRES);
812:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
813:                                     "KSPGMRESSetHapTol_GMRES",
814:                                      KSPGMRESSetHapTol_GMRES);
815:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
816:                                     "KSPGMRESSetCGSRefinementType_GMRES",
817:                                      KSPGMRESSetCGSRefinementType_GMRES);

819:   gmres->haptol              = 1.0e-30;
820:   gmres->q_preallocate       = 0;
821:   gmres->delta_allocate      = GMRES_DELTA_DIRECTIONS;
822:   gmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
823:   gmres->nrs                 = 0;
824:   gmres->sol_temp            = 0;
825:   gmres->max_k               = GMRES_DEFAULT_MAXK;
826:   gmres->Rsvd                = 0;
827:   gmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
828:   return(0);
829: }
830: EXTERN_C_END