Actual source code: gmres.c

petsc-3.4.5 2014-06-29
  2: /*
  3:     This file implements GMRES (a Generalized Minimal Residual) method.
  4:     Reference:  Saad and Schultz, 1986.


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

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

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

 31: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>       /*I  "petscksp.h"  I*/
 32: #define GMRES_DELTA_DIRECTIONS 10
 33: #define GMRES_DEFAULT_MAXK     30
 34: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP,PetscInt,PetscBool,PetscReal*);
 35: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 39: PetscErrorCode    KSPSetUp_GMRES(KSP ksp)
 40: {
 41:   PetscInt       hh,hes,rs,cc;
 43:   PetscInt       max_k,k;
 44:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

 47:   max_k = gmres->max_k;          /* restart size */
 48:   hh    = (max_k + 2) * (max_k + 1);
 49:   hes   = (max_k + 1) * (max_k + 1);
 50:   rs    = (max_k + 2);
 51:   cc    = (max_k + 1);

 53:   PetscMalloc5(hh,PetscScalar,&gmres->hh_origin,hes,PetscScalar,&gmres->hes_origin,rs,PetscScalar,&gmres->rs_origin,cc,PetscScalar,&gmres->cc_origin,cc,PetscScalar,&gmres->ss_origin);
 54:   PetscMemzero(gmres->hh_origin,hh*sizeof(PetscScalar));
 55:   PetscMemzero(gmres->hes_origin,hes*sizeof(PetscScalar));
 56:   PetscMemzero(gmres->rs_origin,rs*sizeof(PetscScalar));
 57:   PetscMemzero(gmres->cc_origin,cc*sizeof(PetscScalar));
 58:   PetscMemzero(gmres->ss_origin,cc*sizeof(PetscScalar));
 59:   PetscLogObjectMemory(ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));

 61:   if (ksp->calc_sings) {
 62:     /* Allocate workspace to hold Hessenberg matrix needed by lapack */
 63:     PetscMalloc((max_k + 3)*(max_k + 9)*sizeof(PetscScalar),&gmres->Rsvd);
 64:     PetscLogObjectMemory(ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
 65:     PetscMalloc(6*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
 66:     PetscLogObjectMemory(ksp,6*(max_k+2)*sizeof(PetscReal));
 67:   }

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

 73:   PetscMalloc((gmres->vecs_allocated)*sizeof(Vec),&gmres->vecs);
 74:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(Vec*),&gmres->user_work);
 75:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&gmres->mwork_alloc);
 76:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(sizeof(Vec*)+sizeof(PetscInt)) + gmres->vecs_allocated*sizeof(Vec));

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

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

 84:     gmres->mwork_alloc[0] = gmres->vv_allocated;
 85:     gmres->nwork_alloc    = 1;
 86:     for (k=0; k<gmres->vv_allocated; k++) {
 87:       gmres->vecs[k] = gmres->user_work[0][k];
 88:     }
 89:   } else {
 90:     gmres->vv_allocated = 5;

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

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:  */
123: PetscErrorCode KSPGMRESCycle(PetscInt *itcount,KSP ksp)
124: {
125:   KSP_GMRES      *gmres = (KSP_GMRES*)(ksp->data);
126:   PetscReal      res_norm,res,hapbnd,tt;
128:   PetscInt       it     = 0, max_k = gmres->max_k;
129:   PetscBool      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:   PetscObjectAMSTakeAccess((PetscObject)ksp);
138:   ksp->rnorm = res;
139:   PetscObjectAMSGrantAccess((PetscObject)ksp);
140:   gmres->it  = (it - 1);
141:   KSPLogResidualHistory(ksp,res);
142:   KSPMonitor(ksp,ksp->its,res);
143:   if (!res) {
144:     if (itcount) *itcount = 0;
145:     ksp->reason = KSP_CONVERGED_ATOL;
146:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
147:     return(0);
148:   }

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

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

165:     /* vv(i+1) . vv(i+1) */
166:     VecNormalize(VEC_VV(it+1),&tt);

168:     /* save the magnitude */
169:     *HH(it+1,it)  = tt;
170:     *HES(it+1,it) = tt;

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

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

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:         if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %G",res);
193:         else {
194:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
195:           break;
196:         }
197:       }
198:     }
199:   }

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

207:   if (itcount) *itcount = it;


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

222: PetscErrorCode KSPSolve_GMRES(KSP ksp)
223: {
225:   PetscInt       its,itcount;
226:   KSP_GMRES      *gmres     = (KSP_GMRES*)ksp->data;
227:   PetscBool      guess_zero = ksp->guess_zero;

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

232:   PetscObjectAMSTakeAccess((PetscObject)ksp);
233:   ksp->its = 0;
234:   PetscObjectAMSGrantAccess((PetscObject)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:     KSPGMRESCycle(&its,ksp);
241:     itcount += its;
242:     if (itcount >= ksp->max_it) {
243:       if (!ksp->reason) 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 KSPReset_GMRES(KSP ksp)
255: {
256:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
258:   PetscInt       i;

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

264:   /* free work vectors */
265:   PetscFree(gmres->vecs);
266:   for (i=0; i<gmres->nwork_alloc; i++) {
267:     VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
268:   }
269:   gmres->nwork_alloc = 0;

271:   PetscFree(gmres->user_work);
272:   PetscFree(gmres->mwork_alloc);
273:   PetscFree(gmres->nrs);
274:   VecDestroy(&gmres->sol_temp);
275:   PetscFree(gmres->Rsvd);
276:   PetscFree(gmres->Dsvd);
277:   PetscFree(gmres->orthogwork);

279:   gmres->sol_temp       = 0;
280:   gmres->vv_allocated   = 0;
281:   gmres->vecs_allocated = 0;
282:   gmres->sol_temp       = 0;
283:   return(0);
284: }

288: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
289: {

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

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

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

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

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

340:     PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D GRS(it) = %G",it,PetscAbsScalar(*GRS(it)));
341:     return(0);
342:   }
343:   for (ii=1; ii<=it; ii++) {
344:     k  = it - ii;
345:     tt = *GRS(k);
346:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
347:     if (*HH(k,k) == 0.0) {
348:       ksp->reason = KSP_DIVERGED_BREAKDOWN;

350:       PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D",k);
351:       return(0);
352:     }
353:     nrs[k] = tt / *HH(k,k);
354:   }

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

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

380:   hh = HH(0,it);
381:   cc = CC(0);
382:   ss = SS(0);

384:   /* Apply all the previously computed plane rotations to the new column
385:      of the Hessenberg matrix */
386:   for (j=1; j<=it; j++) {
387:     tt  = *hh;
388:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
389:     hh++;
390:     *hh = *cc++ * *hh - (*ss++ * tt);
391:   }

393:   /*
394:     compute the new plane rotation, and apply it to:
395:      1) the right-hand-side of the Hessenberg system
396:      2) the new column of the Hessenberg matrix
397:     thus obtaining the updated value of the residual
398:   */
399:   if (!hapend) {
400:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
401:     if (tt == 0.0) {
402:       ksp->reason = KSP_DIVERGED_NULL;
403:       return(0);
404:     }
405:     *cc        = *hh / tt;
406:     *ss        = *(hh+1) / tt;
407:     *GRS(it+1) = -(*ss * *GRS(it));
408:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
409:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
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". */

419:     *res = 0.0;
420:   }
421:   return(0);
422: }
423: /*
424:    This routine allocates more work vectors, starting from VEC_VV(it).
425:  */
428: PetscErrorCode KSPGMRESGetNewVectors(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;

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

448:   gmres->mwork_alloc[nwork] = nalloc;
449:   for (k=0; k<nalloc; k++) {
450:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
451:   }
452:   gmres->nwork_alloc++;
453:   return(0);
454: }

458: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
459: {
460:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

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

477:   KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
478:   if (result) *result = ptr;
479:   return(0);
480: }

484: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
485: {
486:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
487:   const char     *cstr;
489:   PetscBool      iascii,isstring;

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

524: /*@C
525:    KSPGMRESMonitorKrylov - 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: KSPMonitorSet(), KSPMonitorDefault(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
541: @*/
542: PetscErrorCode  KSPGMRESMonitorKrylov(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;
549:   PetscBool      flg;

552:   PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
553:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
554:   if (!flg) {
555:     PetscViewerSetType(viewer,PETSCVIEWERDRAW);
556:     PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
557:   }

559:   x    = VEC_VV(gmres->it+1);
560:   VecView(x,viewer);
561:   return(0);
562: }

566: PetscErrorCode KSPSetFromOptions_GMRES(KSP ksp)
567: {
569:   PetscInt       restart;
570:   PetscReal      haptol;
571:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
572:   PetscBool      flg;

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

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

610:   if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
611:   gmres->haptol = tol;
612:   return(0);
613: }

617: PetscErrorCode  KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
618: {
619:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

622:   *max_k = gmres->max_k;
623:   return(0);
624: }

628: PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
629: {
630:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

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

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

657: PetscErrorCode  KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
658: {
660:   *fcn = ((KSP_GMRES*)ksp->data)->orthog;
661:   return(0);
662: }

666: PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
667: {
668:   KSP_GMRES *gmres;

671:   gmres = (KSP_GMRES*)ksp->data;
672:   gmres->q_preallocate = 1;
673:   return(0);
674: }

678: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
679: {
680:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

683:   gmres->cgstype = type;
684:   return(0);
685: }

689: PetscErrorCode  KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
690: {
691:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

694:   *type = gmres->cgstype;
695:   return(0);
696: }

700: /*@
701:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
702:          in the classical Gram Schmidt orthogonalization.

704:    Logically Collective on KSP

706:    Input Parameters:
707: +  ksp - the Krylov space context
708: -  type - the type of refinement

710:   Options Database:
711: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

713:    Level: intermediate

715: .keywords: KSP, GMRES, iterative refinement

717: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
718:           KSPGMRESGetOrthogonalization()
719: @*/
720: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
721: {

727:   PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
728:   return(0);
729: }

733: /*@
734:    KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
735:          in the classical Gram Schmidt orthogonalization.

737:    Not Collective

739:    Input Parameter:
740: .  ksp - the Krylov space context

742:    Output Parameter:
743: .  type - the type of refinement

745:   Options Database:
746: .  -ksp_gmres_cgs_refinement_type <never,ifneeded,always>

748:    Level: intermediate

750: .keywords: KSP, GMRES, iterative refinement

752: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
753:           KSPGMRESGetOrthogonalization()
754: @*/
755: PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
756: {

761:   PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
762:   return(0);
763: }


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

771:    Logically Collective on KSP

773:    Input Parameters:
774: +  ksp - the Krylov space context
775: -  restart - integer restart value

777:   Options Database:
778: .  -ksp_gmres_restart <positive integer>

780:     Note: The default value is 30.

782:    Level: intermediate

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

786: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
787: @*/
788: PetscErrorCode  KSPGMRESSetRestart(KSP ksp, PetscInt restart)
789: {


795:   PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
796:   return(0);
797: }

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

804:    Not Collective

806:    Input Parameter:
807: .  ksp - the Krylov space context

809:    Output Parameter:
810: .   restart - integer restart value

812:     Note: The default value is 30.

814:    Level: intermediate

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

818: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
819: @*/
820: PetscErrorCode  KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
821: {

825:   PetscTryMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
826:   return(0);
827: }

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

834:    Logically Collective on KSP

836:    Input Parameters:
837: +  ksp - the Krylov space context
838: -  tol - the tolerance

840:   Options Database:
841: .  -ksp_gmres_haptol <positive real value>

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

847:    Level: intermediate

849: .keywords: KSP, GMRES, tolerance

851: .seealso: KSPSetTolerances()
852: @*/
853: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
854: {

859:   PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
860:   return(0);
861: }

863: /*MC
864:      KSPGMRES - Implements the Generalized Minimal Residual method.
865:                 (Saad and Schultz, 1986) with restart


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

879:    Level: beginner

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

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

887: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
888:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
889:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
890:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()

892: M*/

896: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
897: {
898:   KSP_GMRES      *gmres;

902:   PetscNewLog(ksp,KSP_GMRES,&gmres);
903:   ksp->data = (void*)gmres;

905:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
906:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);

908:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;
909:   ksp->ops->setup                        = KSPSetUp_GMRES;
910:   ksp->ops->solve                        = KSPSolve_GMRES;
911:   ksp->ops->reset                        = KSPReset_GMRES;
912:   ksp->ops->destroy                      = KSPDestroy_GMRES;
913:   ksp->ops->view                         = KSPView_GMRES;
914:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
915:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
916:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

918:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
919:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
920:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
921:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
922:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
923:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
924:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
925:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

927:   gmres->haptol         = 1.0e-30;
928:   gmres->q_preallocate  = 0;
929:   gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
930:   gmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
931:   gmres->nrs            = 0;
932:   gmres->sol_temp       = 0;
933:   gmres->max_k          = GMRES_DEFAULT_MAXK;
934:   gmres->Rsvd           = 0;
935:   gmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
936:   gmres->orthogwork     = 0;
937:   return(0);
938: }