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/gmresimpl.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       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);

 59:   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);
 60:   PetscMemzero(gmres->hh_origin,hh*sizeof(PetscScalar));
 61:   PetscMemzero(gmres->hes_origin,hes*sizeof(PetscScalar));
 62:   PetscMemzero(gmres->rs_origin,rs*sizeof(PetscScalar));
 63:   PetscMemzero(gmres->cc_origin,cc*sizeof(PetscScalar));
 64:   PetscMemzero(gmres->ss_origin,cc*sizeof(PetscScalar));
 65:   PetscLogObjectMemory(ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));

 67:   if (ksp->calc_sings) {
 68:     /* Allocate workspace to hold Hessenberg matrix needed by lapack */
 69:     PetscMalloc((max_k + 3)*(max_k + 9)*sizeof(PetscScalar),&gmres->Rsvd);
 70:     PetscLogObjectMemory(ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
 71:     PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
 72:     PetscLogObjectMemory(ksp,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(PetscInt),&gmres->mwork_alloc);
 81:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));

 83:   if (gmres->q_preallocate) {
 84:     gmres->vv_allocated   = VEC_OFFSET + 2 + max_k;
 85:     KSPGetVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,PETSC_NULL);
 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:     KSPGetVecs(ksp,5,&gmres->user_work[0],0,PETSC_NULL);
 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: PetscErrorCode GMREScycle(PetscInt *itcount,KSP ksp)
125: {
126:   KSP_GMRES      *gmres = (KSP_GMRES *)(ksp->data);
127:   PetscReal      res_norm,res,hapbnd,tt;
129:   PetscInt       it = 0, max_k = gmres->max_k;
130:   PetscTruth     hapend = PETSC_FALSE;

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

137:   /* check for the convergence */
138:   PetscObjectTakeAccess(ksp);
139:   ksp->rnorm = res;
140:   PetscObjectGrantAccess(ksp);
141:   gmres->it = (it - 1);
142:   KSPLogResidualHistory(ksp,res);
143:   KSPMonitor(ksp,ksp->its,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:     if (it) {
154:       KSPLogResidualHistory(ksp,res);
155:       KSPMonitor(ksp,ksp->its,res);
156:     }
157:     gmres->it = (it - 1);
158:     if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
159:       GMRESGetNewVectors(ksp,it+1);
160:     }
161:     KSP_PCApplyBAorAB(ksp,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);

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

166:     /* vv(i+1) . vv(i+1) */
167:     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 = %G tt = %G\n",hapbnd,tt);
177:       hapend = PETSC_TRUE;
178:     }
179:     GMRESUpdateHessenberg(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:         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 (it && (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 && ksp->pc_side != PC_RIGHT) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Use right preconditioning -ksp_right_pc if want unpreconditioned norm)");

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:       if (!ksp->reason) 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 matrices */
263:   PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_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:     VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
332:     return(0);
333:   }
334:   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)));
335:   if (*HH(it,it) != 0.0) {
336:     nrs[it] = *GRS(it) / *HH(it,it);
337:   } else {
338:     nrs[it] = 0.0;
339:   }
340:   for (ii=1; ii<=it; ii++) {
341:     k   = it - ii;
342:     tt  = *GRS(k);
343:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
344:     if (*HH(k,k) == 0.0) SETERRQ2(PETSC_ERR_CONV_FAILED,"HH(k,k) is identically zero; it = %D k = %D",it,k);
345:     nrs[k]   = tt / *HH(k,k);
346:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

533:    Collective on KSP

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

541:    Level: intermediate

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

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

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

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

562:   return(0);
563: }

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

576:   PetscOptionsHead("KSP GMRES Options");
577:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
578:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
579:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
580:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
581:     flg  = PETSC_FALSE;
582:     PetscOptionsTruth("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,PETSC_NULL);
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:     flg  = PETSC_FALSE;
591:     PetscOptionsTruth("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,PETSC_NULL);
592:     if (flg) {
593:       PetscViewers viewers;
594:       PetscViewersCreate(((PetscObject)ksp)->comm,&viewers);
595:       KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void*))PetscViewersDestroy);
596:     }
597:   PetscOptionsTail();
598:   return(0);
599: }

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


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

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

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

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

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

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

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

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

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

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

688:    Collective on KSP

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

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

697:    Level: intermediate

699: .keywords: KSP, GMRES, iterative refinement

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

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

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

721:    Collective on KSP

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

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

730:     Note: The default value is 30.

732:    Level: intermediate

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

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

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

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

752:    Collective on KSP

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

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

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

765:    Level: intermediate

767: .keywords: KSP, GMRES, tolerance

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

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

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


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

796:    Level: beginner

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

802: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
803:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
804:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
805:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov()

807: M*/

812: PetscErrorCode  KSPCreate_GMRES(KSP ksp)
813: {
814:   KSP_GMRES      *gmres;

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


822:   ksp->normtype                          = KSP_NORM_PRECONDITIONED;
823:   ksp->pc_side                           = PC_LEFT;

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

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

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