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/sles/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);

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

 46: #if 0
 47:   if (ksp->pc_side == PC_SYMMETRIC) {
 48:     SETERRQ(2,"no symmetric preconditioning for KSPGMRES");
 49:   }
 50: #endif
 51:   max_k         = gmres->max_k;
 52:   hh            = (max_k + 2) * (max_k + 1);
 53:   hes           = (max_k + 1) * (max_k + 1);
 54:   rs            = (max_k + 2);
 55:   cc            = (max_k + 1);
 56:   size          = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);

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

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

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

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

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

108: .        gmres  - structure containing parameters and work areas

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

130:   ierr    = VecNorm(VEC_VV(0),NORM_2,&res_norm);
131:   res     = res_norm;
132:   *GRS(0) = res_norm;

134:   /* check for the convergence */
135:   if (!res) {
136:     if (itcount) *itcount = 0;
137:     ksp->reason = KSP_CONVERGED_ATOL;
138:     PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entryn");
139:     return(0);
140:   }

142:   /* scale VEC_VV (the initial residual) */
143:   tmp = 1.0/res_norm; VecScale(&tmp,VEC_VV(0));

145:   PetscObjectTakeAccess(ksp);
146:   ksp->rnorm = res;
147:   PetscObjectGrantAccess(ksp);
148:   gmres->it = (it - 1);
149:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
150:   while (!ksp->reason && it < max_k && ksp->its < 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:     VecNorm(VEC_VV(it+1),NORM_2,&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:       tmp = 1.0/tt; VecScale(&tmp,VEC_VV(it+1));
173:     } else {
174:       PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %gn",hapbnd,tt);
175:       hapend = PETSC_TRUE;
176:     }
177:     GMRESUpdateHessenberg(ksp,it,hapend,&res);
178:     it++;
179:     gmres->it  = (it-1);  /* For converged */
180:     PetscObjectTakeAccess(ksp);
181:     ksp->its++;
182:     ksp->rnorm = res;
183:     PetscObjectGrantAccess(ksp);

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

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

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

203:   if (itcount) *itcount    = it;


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

214:   return(0);
215: }

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

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

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

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

249: int KSPDestroy_GMRES(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);}
271:   PetscFree(gmres);

273:   return(0);
274: }
275: /*
276:     BuildGmresSoln - create the solution from the starting vector and the
277:     current iterates.

279:     Input parameters:
280:         nrs - work area of size it + 1.
281:         vs  - index of initial guess
282:         vdest - index of result.  Note that vs may == vdest (replace
283:                 guess with the solution).

285:      This is an internal routine that knows about the GMRES internals.
286:  */
287: static int BuildGmresSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,int it)
288: {
289:   PetscScalar tt,zero = 0.0,one = 1.0;
290:   int         ierr,ii,k,j;
291:   KSP_GMRES   *gmres = (KSP_GMRES *)(ksp->data);

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

296:   /* If it is < 0, no gmres steps have been performed */
297:   if (it < 0) {
298:     if (vdest != vs) {
299:       VecCopy(vs,vdest);
300:     }
301:     return(0);
302:   }
303:   if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
304:   if (*HH(it,it) != 0.0) {
305:     nrs[it] = *GRS(it) / *HH(it,it);
306:   } else {
307:     nrs[it] = 0.0;
308:   }
309:   for (ii=1; ii<=it; ii++) {
310:     k   = it - ii;
311:     tt  = *GRS(k);
312:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
313:     nrs[k]   = tt / *HH(k,k);
314:   }

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

320:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
321:   /* add solution to previous solution */
322:   if (vdest != vs) {
323:     VecCopy(vs,vdest);
324:   }
325:   VecAXPY(&one,VEC_TEMP,vdest);
326:   return(0);
327: }
328: /*
329:    Do the scalar work for the orthogonalization.  Return new residual.
330:  */
331: static int GMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
332: {
333:   PetscScalar *hh,*cc,*ss,tt;
334:   int         j;
335:   KSP_GMRES   *gmres = (KSP_GMRES *)(ksp->data);

338:   hh  = HH(0,it);
339:   cc  = CC(0);
340:   ss  = SS(0);

342:   /* Apply all the previously computed plane rotations to the new column
343:      of the Hessenberg matrix */
344:   for (j=1; j<=it; j++) {
345:     tt  = *hh;
346: #if defined(PETSC_USE_COMPLEX)
347:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
348: #else
349:     *hh = *cc * tt + *ss * *(hh+1);
350: #endif
351:     hh++;
352:     *hh = *cc++ * *hh - (*ss++ * tt);
353:   }

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

400:   nalloc = gmres->delta_allocate;
401:   /* Adjust the number to allocate to make sure that we don't exceed the
402:     number of available slots */
403:   if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
404:     nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
405:   }
406:   if (!nalloc) return(0);

408:   gmres->vv_allocated += nalloc;
409:   VecDuplicateVecs(ksp->vec_rhs,nalloc,&gmres->user_work[nwork]);
410:   PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
411:   gmres->mwork_alloc[nwork] = nalloc;
412:   for (k=0; k<nalloc; k++) {
413:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
414:   }
415:   gmres->nwork_alloc++;
416:   return(0);
417: }

419: int KSPBuildSolution_GMRES(KSP ksp,Vec  ptr,Vec *result)
420: {
421:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
422:   int       ierr;

425:   if (!ptr) {
426:     if (!gmres->sol_temp) {
427:       VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
428:       PetscLogObjectParent(ksp,gmres->sol_temp);
429:     }
430:     ptr = gmres->sol_temp;
431:   }
432:   if (!gmres->nrs) {
433:     /* allocate the work area */
434:     PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
435:     PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
436:   }

438:   BuildGmresSoln(gmres->nrs,VEC_SOLN,ptr,ksp,gmres->it);
439:   *result = ptr;
440:   return(0);
441: }

443: int KSPView_GMRES(KSP ksp,PetscViewer viewer)
444: {
445:   KSP_GMRES  *gmres = (KSP_GMRES *)ksp->data;
446:   char       *cstr;
447:   int        ierr;
448:   PetscTruth isascii,isstring;

451:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
452:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
453:   if (gmres->orthog == KSPGMRESUnmodifiedGramSchmidtOrthogonalization) {
454:     cstr = "Unmodified Gram-Schmidt Orthogonalization";
455:   } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
456:     cstr = "Modified Gram-Schmidt Orthogonalization";
457:   } else if (gmres->orthog == KSPGMRESIROrthogonalization) {
458:     cstr = "Unmodified Gram-Schmidt + 1 step Iterative Refinement Orthogonalization";
459:   } else {
460:     cstr = "unknown orthogonalization";
461:   }
462:   if (isascii) {
463:     PetscViewerASCIIPrintf(viewer,"  GMRES: restart=%d, using %sn",gmres->max_k,cstr);
464:     PetscViewerASCIIPrintf(viewer,"  GMRES: happy breakdown tolerance %gn",gmres->haptol);
465:   } else if (isstring) {
466:     PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,gmres->max_k);
467:   } else {
468:     SETERRQ1(1,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
469:   }
470:   return(0);
471: }

473: /*@C
474:    KSPGMRESKrylovMonitor - Calls VecView() for each direction in the 
475:    GMRES accumulated Krylov space.

477:    Collective on KSP

479:    Input Parameters:
480: +  ksp - the KSP context
481: .  its - iteration number
482: .  fgnorm - 2-norm of residual (or gradient)
483: -  a viewers object created with PetscViewersCreate()

485:    Level: intermediate

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

489: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
490: @*/
491: int KSPGMRESKrylovMonitor(KSP ksp,int its,PetscReal fgnorm,void *dummy)
492: {
493:   PetscViewers viewers = (PetscViewers)dummy;
494:   KSP_GMRES    *gmres = (KSP_GMRES*)ksp->data;
495:   int          ierr;
496:   Vec          x;
497:   PetscViewer  viewer;

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

503:   x      = VEC_VV(gmres->it+1);
504:   ierr   = VecView(x,viewer);

506:   return(0);
507: }

509: int KSPSetFromOptions_GMRES(KSP ksp)
510: {
511:   int        ierr,restart;
512:   PetscReal  haptol;
513:   KSP_GMRES  *gmres = (KSP_GMRES*)ksp->data;
514:   PetscTruth flg;

517:   PetscOptionsHead("KSP GMRES Options");
518:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
519:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
520:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
521:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
522:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
523:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
524:     PetscOptionsLogicalGroupBegin("-ksp_gmres_unmodifiedgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
525:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESUnmodifiedGramSchmidtOrthogonalization);}
526:     PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
527:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
528:     PetscOptionsLogicalGroupEnd("-ksp_gmres_irorthog","Classical Gram-Schmidt + iterative refinement","KSPGMRESSetOrthogonalization",&flg);
529:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESIROrthogonalization);}
530:     PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
531:     if (flg) {
532:       PetscViewers viewers;
533:       PetscViewersCreate(ksp->comm,&viewers);
534:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
535:     }
536:   PetscOptionsTail();
537:   return(0);
538: }

540: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
541: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);


544: EXTERN_C_BEGIN
545: int KSPGMRESSetHapTol_GMRES(KSP ksp,double tol)
546: {
547:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;

550:   if (tol < 0.0) SETERRQ(1,"Tolerance must be non-negative");
551:   gmres->haptol = tol;
552:   return(0);
553: }
554: EXTERN_C_END

556: EXTERN_C_BEGIN
557: int KSPGMRESSetRestart_GMRES(KSP ksp,int max_k)
558: {
559:   KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;

562:   if (max_k < 1) SETERRQ(1,"Restart must be positive");
563:   gmres->max_k = max_k;
564:   return(0);
565: }
566: EXTERN_C_END

568: EXTERN_C_BEGIN
569: int KSPGMRESSetOrthogonalization_GMRES(KSP ksp,int (*fcn)(KSP,int))
570: {
573:   ((KSP_GMRES *)ksp->data)->orthog = fcn;
574:   return(0);
575: }
576: EXTERN_C_END

578: EXTERN_C_BEGIN
579: int KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
580: {
581:   KSP_GMRES *gmres;

584:   gmres = (KSP_GMRES *)ksp->data;
585:   gmres->q_preallocate = 1;
586:   return(0);
587: }
588: EXTERN_C_END

590: EXTERN_C_BEGIN
591: int KSPCreate_GMRES(KSP ksp)
592: {
593:   KSP_GMRES *gmres;
594:   int       ierr;

597:   PetscNew(KSP_GMRES,&gmres);
598:   ierr  = PetscMemzero(gmres,sizeof(KSP_GMRES));
599:   PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
600:   ksp->data                              = (void*)gmres;
601:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;

603:   ksp->ops->setup                        = KSPSetUp_GMRES;
604:   ksp->ops->solve                        = KSPSolve_GMRES;
605:   ksp->ops->destroy                      = KSPDestroy_GMRES;
606:   ksp->ops->view                         = KSPView_GMRES;
607:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
608:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
609:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

611:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
612:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
613:                                      KSPGMRESSetPreAllocateVectors_GMRES);
614:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
615:                                     "KSPGMRESSetOrthogonalization_GMRES",
616:                                      KSPGMRESSetOrthogonalization_GMRES);
617:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
618:                                     "KSPGMRESSetRestart_GMRES",
619:                                      KSPGMRESSetRestart_GMRES);
620:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
621:                                     "KSPGMRESSetHapTol_GMRES",
622:                                      KSPGMRESSetHapTol_GMRES);

624:   gmres->haptol              = 1.0e-30;
625:   gmres->q_preallocate       = 0;
626:   gmres->delta_allocate      = GMRES_DELTA_DIRECTIONS;
627:   gmres->orthog              = KSPGMRESUnmodifiedGramSchmidtOrthogonalization;
628:   gmres->nrs                 = 0;
629:   gmres->sol_temp            = 0;
630:   gmres->max_k               = GMRES_DEFAULT_MAXK;
631:   gmres->Rsvd                = 0;
632:   return(0);
633: }
634: EXTERN_C_END