Actual source code: gmres.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

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

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

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

 51:   PetscCalloc5(hh,&gmres->hh_origin,hes,&gmres->hes_origin,rs,&gmres->rs_origin,cc,&gmres->cc_origin,cc,&gmres->ss_origin);
 52:   PetscLogObjectMemory((PetscObject)ksp,(hh + hes + rs + 2*cc)*sizeof(PetscScalar));

 54:   if (ksp->calc_sings) {
 55:     /* Allocate workspace to hold Hessenberg matrix needed by lapack */
 56:     PetscMalloc1((max_k + 3)*(max_k + 9),&gmres->Rsvd);
 57:     PetscLogObjectMemory((PetscObject)ksp,(max_k + 3)*(max_k + 9)*sizeof(PetscScalar));
 58:     PetscMalloc1(6*(max_k+2),&gmres->Dsvd);
 59:     PetscLogObjectMemory((PetscObject)ksp,6*(max_k+2)*sizeof(PetscReal));
 60:   }

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

 66:   PetscMalloc1(gmres->vecs_allocated,&gmres->vecs);
 67:   PetscMalloc1(VEC_OFFSET+2+max_k,&gmres->user_work);
 68:   PetscMalloc1(VEC_OFFSET+2+max_k,&gmres->mwork_alloc);
 69:   PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+2+max_k)*(sizeof(Vec*)+sizeof(PetscInt)) + gmres->vecs_allocated*sizeof(Vec));

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

 74:     KSPCreateVecs(ksp,gmres->vv_allocated,&gmres->user_work[0],0,NULL);
 75:     PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);

 77:     gmres->mwork_alloc[0] = gmres->vv_allocated;
 78:     gmres->nwork_alloc    = 1;
 79:     for (k=0; k<gmres->vv_allocated; k++) {
 80:       gmres->vecs[k] = gmres->user_work[0][k];
 81:     }
 82:   } else {
 83:     gmres->vv_allocated = 5;

 85:     KSPCreateVecs(ksp,5,&gmres->user_work[0],0,NULL);
 86:     PetscLogObjectParents(ksp,5,gmres->user_work[0]);

 88:     gmres->mwork_alloc[0] = 5;
 89:     gmres->nwork_alloc    = 1;
 90:     for (k=0; k<gmres->vv_allocated; k++) {
 91:       gmres->vecs[k] = gmres->user_work[0][k];
 92:     }
 93:   }
 94:   return(0);
 95: }

 97: /*
 98:     Run gmres, possibly with restart.  Return residual history if requested.
 99:     input parameters:

101: .        gmres  - structure containing parameters and work areas

103:     output parameters:
104: .        nres    - residuals (from preconditioned system) at each step.
105:                   If restarting, consider passing nres+it.  If null,
106:                   ignored
107: .        itcount - number of iterations used.  nres[0] to nres[itcount]
108:                   are defined.  If null, ignored.

110:     Notes:
111:     On entry, the value in vector VEC_VV(0) should be the initial residual
112:     (this allows shortcuts where the initial preconditioned residual is 0).
113:  */
114: PetscErrorCode KSPGMRESCycle(PetscInt *itcount,KSP ksp)
115: {
116:   KSP_GMRES      *gmres = (KSP_GMRES*)(ksp->data);
117:   PetscReal      res,hapbnd,tt;
119:   PetscInt       it     = 0, max_k = gmres->max_k;
120:   PetscBool      hapend = PETSC_FALSE;

123:   if (itcount) *itcount = 0;
124:   VecNormalize(VEC_VV(0),&res);
125:   KSPCheckNorm(ksp,res);

127:   /* the constant .1 is arbitrary, just some measure at how incorrect the residuals are */
128:   if ((ksp->rnorm > 0.0) && (PetscAbsReal(res-ksp->rnorm) > gmres->breakdowntol*gmres->rnorm0)) {
129:       if (ksp->errorifnotconverged) SETERRQ3(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Residual norm computed by GMRES recursion formula %g is far from the computed residual norm %g at restart, residual norm at start of cycle %g",(double)ksp->rnorm,(double)res,(double)gmres->rnorm0);
130:       else {
131:         PetscInfo3(ksp,"Residual norm computed by GMRES recursion formula %g is far from the computed residual norm %g at restart, residual norm at start of cycle %g",(double)ksp->rnorm,(double)res,(double)gmres->rnorm0);
132:         ksp->reason =  KSP_DIVERGED_BREAKDOWN;
133:         return(0);
134:       }
135:   }
136:   *GRS(0) = gmres->rnorm0 = res;

138:   /* check for the convergence */
139:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
140:   ksp->rnorm = res;
141:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
142:   gmres->it  = (it - 1);
143:   KSPLogResidualHistory(ksp,res);
144:   KSPMonitor(ksp,ksp->its,res);
145:   if (!res) {
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:       KSPGMRESGetNewVectors(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);
165:     if (ksp->reason) break;

167:     /* vv(i+1) . vv(i+1) */
168:     VecNormalize(VEC_VV(it+1),&tt);
169:     KSPCheckNorm(ksp,tt);

171:     /* save the magnitude */
172:     *HH(it+1,it)  = tt;
173:     *HES(it+1,it) = tt;

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

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

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

192:     /* Catch error in happy breakdown and signal convergence and break from loop */
193:     if (hapend) {
194:       if (ksp->normtype == KSP_NORM_NONE) { /* convergence test was skipped in this case */
195:         ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
196:       } else if (!ksp->reason) {
197:         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",(double)res);
198:         else {
199:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
200:           break;
201:         }
202:       }
203:     }
204:   }

206:   /* Monitor if we know that we will not return for a restart */
207:   if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
208:     KSPLogResidualHistory(ksp,res);
209:     KSPMonitor(ksp,ksp->its,res);
210:   }

212:   if (itcount) *itcount = it;


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

225: PetscErrorCode KSPSolve_GMRES(KSP ksp)
226: {
228:   PetscInt       its,itcount,i;
229:   KSP_GMRES      *gmres     = (KSP_GMRES*)ksp->data;
230:   PetscBool      guess_zero = ksp->guess_zero;
231:   PetscInt       N = gmres->max_k + 1;

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

236:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
237:   ksp->its = 0;
238:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

240:   itcount     = 0;
241:   gmres->fullcycle = 0;
242:   ksp->reason = KSP_CONVERGED_ITERATING;
243:   ksp->rnorm  = -1.0; /* special marker for KSPGMRESCycle() */
244:   while (!ksp->reason) {
245:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
246:     KSPGMRESCycle(&its,ksp);
247:     /* Store the Hessenberg matrix and the basis vectors of the Krylov subspace
248:     if the cycle is complete for the computation of the Ritz pairs */
249:     if (its == gmres->max_k) {
250:       gmres->fullcycle++;
251:       if (ksp->calc_ritz) {
252:         if (!gmres->hes_ritz) {
253:           PetscMalloc1(N*N,&gmres->hes_ritz);
254:           PetscLogObjectMemory((PetscObject)ksp,N*N*sizeof(PetscScalar));
255:           VecDuplicateVecs(VEC_VV(0),N,&gmres->vecb);
256:         }
257:         PetscArraycpy(gmres->hes_ritz,gmres->hes_origin,N*N);
258:         for (i=0; i<gmres->max_k+1; i++) {
259:           VecCopy(VEC_VV(i),gmres->vecb[i]);
260:         }
261:       }
262:     }
263:     itcount += its;
264:     if (itcount >= ksp->max_it) {
265:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
266:       break;
267:     }
268:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
269:   }
270:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
271:   return(0);
272: }

274: PetscErrorCode KSPReset_GMRES(KSP ksp)
275: {
276:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
278:   PetscInt       i;

281:   /* Free the Hessenberg matrices */
282:   PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin);
283:   PetscFree(gmres->hes_ritz);

285:   /* free work vectors */
286:   PetscFree(gmres->vecs);
287:   for (i=0; i<gmres->nwork_alloc; i++) {
288:     VecDestroyVecs(gmres->mwork_alloc[i],&gmres->user_work[i]);
289:   }
290:   gmres->nwork_alloc = 0;
291:   if (gmres->vecb)  {
292:     VecDestroyVecs(gmres->max_k+1,&gmres->vecb);
293:   }

295:   PetscFree(gmres->user_work);
296:   PetscFree(gmres->mwork_alloc);
297:   PetscFree(gmres->nrs);
298:   VecDestroy(&gmres->sol_temp);
299:   PetscFree(gmres->Rsvd);
300:   PetscFree(gmres->Dsvd);
301:   PetscFree(gmres->orthogwork);

303:   gmres->vv_allocated   = 0;
304:   gmres->vecs_allocated = 0;
305:   gmres->sol_temp       = NULL;
306:   return(0);
307: }

309: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
310: {

314:   KSPReset_GMRES(ksp);
315:   PetscFree(ksp->data);
316:   /* clear composed functions */
317:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",NULL);
318:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",NULL);
319:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",NULL);
320:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",NULL);
321:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",NULL);
322:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",NULL);
323:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetBreakdownTolerance_C",NULL);
324:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",NULL);
325:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",NULL);
326:   return(0);
327: }
328: /*
329:     KSPGMRESBuildSoln - create the solution from the starting vector and the
330:     current iterates.

332:     Input parameters:
333:         nrs - work area of size it + 1.
334:         vs  - index of initial guess
335:         vdest - index of result.  Note that vs may == vdest (replace
336:                 guess with the solution).

338:      This is an internal routine that knows about the GMRES internals.
339:  */
340: static PetscErrorCode KSPGMRESBuildSoln(PetscScalar *nrs,Vec vs,Vec vdest,KSP ksp,PetscInt it)
341: {
342:   PetscScalar    tt;
344:   PetscInt       ii,k,j;
345:   KSP_GMRES      *gmres = (KSP_GMRES*)(ksp->data);

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

350:   /* If it is < 0, no gmres steps have been performed */
351:   if (it < 0) {
352:     VecCopy(vs,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
353:     return(0);
354:   }
355:   if (*HH(it,it) != 0.0) {
356:     nrs[it] = *GRS(it) / *HH(it,it);
357:   } else {
358:     if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the break down in GMRES; HH(it,it) = 0");
359:     else ksp->reason = KSP_DIVERGED_BREAKDOWN;

361:     PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D GRS(it) = %g\n",it,(double)PetscAbsScalar(*GRS(it)));
362:     return(0);
363:   }
364:   for (ii=1; ii<=it; ii++) {
365:     k  = it - ii;
366:     tt = *GRS(k);
367:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
368:     if (*HH(k,k) == 0.0) {
369:       if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
370:       else {
371:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
372:         PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
373:         return(0);
374:       }
375:     }
376:     nrs[k] = tt / *HH(k,k);
377:   }

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

383:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
384:   /* add solution to previous solution */
385:   if (vdest != vs) {
386:     VecCopy(vs,vdest);
387:   }
388:   VecAXPY(vdest,1.0,VEC_TEMP);
389:   return(0);
390: }
391: /*
392:    Do the scalar work for the orthogonalization.  Return new residual norm.
393:  */
394: static PetscErrorCode KSPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
395: {
396:   PetscScalar *hh,*cc,*ss,tt;
397:   PetscInt    j;
398:   KSP_GMRES   *gmres = (KSP_GMRES*)(ksp->data);

401:   hh = HH(0,it);
402:   cc = CC(0);
403:   ss = SS(0);

405:   /* Apply all the previously computed plane rotations to the new column
406:      of the Hessenberg matrix */
407:   for (j=1; j<=it; j++) {
408:     tt  = *hh;
409:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
410:     hh++;
411:     *hh = *cc++ * *hh - (*ss++ * tt);
412:   }

414:   /*
415:     compute the new plane rotation, and apply it to:
416:      1) the right-hand-side of the Hessenberg system
417:      2) the new column of the Hessenberg matrix
418:     thus obtaining the updated value of the residual
419:   */
420:   if (!hapend) {
421:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
422:     if (tt == 0.0) {
423:       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"tt == 0.0");
424:       else {
425:         ksp->reason = KSP_DIVERGED_NULL;
426:         return(0);
427:       }
428:     }
429:     *cc        = *hh / tt;
430:     *ss        = *(hh+1) / tt;
431:     *GRS(it+1) = -(*ss * *GRS(it));
432:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
433:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
434:     *res       = PetscAbsScalar(*GRS(it+1));
435:   } else {
436:     /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
437:             another rotation matrix (so RH doesn't change).  The new residual is
438:             always the new sine term times the residual from last time (GRS(it)),
439:             but now the new sine rotation would be zero...so the residual should
440:             be zero...so we will multiply "zero" by the last residual.  This might
441:             not be exactly what we want to do here -could just return "zero". */

443:     *res = 0.0;
444:   }
445:   return(0);
446: }
447: /*
448:    This routine allocates more work vectors, starting from VEC_VV(it).
449:  */
450: PetscErrorCode KSPGMRESGetNewVectors(KSP ksp,PetscInt it)
451: {
452:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
454:   PetscInt       nwork = gmres->nwork_alloc,k,nalloc;

457:   nalloc = PetscMin(ksp->max_it,gmres->delta_allocate);
458:   /* Adjust the number to allocate to make sure that we don't exceed the
459:     number of available slots */
460:   if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated) {
461:     nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
462:   }
463:   if (!nalloc) return(0);

465:   gmres->vv_allocated += nalloc;

467:   KSPCreateVecs(ksp,nalloc,&gmres->user_work[nwork],0,NULL);
468:   PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);

470:   gmres->mwork_alloc[nwork] = nalloc;
471:   for (k=0; k<nalloc; k++) {
472:     gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
473:   }
474:   gmres->nwork_alloc++;
475:   return(0);
476: }

478: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
479: {
480:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

484:   if (!ptr) {
485:     if (!gmres->sol_temp) {
486:       VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
487:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)gmres->sol_temp);
488:     }
489:     ptr = gmres->sol_temp;
490:   }
491:   if (!gmres->nrs) {
492:     /* allocate the work area */
493:     PetscMalloc1(gmres->max_k,&gmres->nrs);
494:     PetscLogObjectMemory((PetscObject)ksp,gmres->max_k);
495:   }

497:   KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
498:   if (result) *result = ptr;
499:   return(0);
500: }

502: PetscErrorCode KSPView_GMRES(KSP ksp,PetscViewer viewer)
503: {
504:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
505:   const char     *cstr;
507:   PetscBool      iascii,isstring;

510:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
511:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
512:   if (gmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
513:     switch (gmres->cgstype) {
514:     case (KSP_GMRES_CGS_REFINE_NEVER):
515:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
516:       break;
517:     case (KSP_GMRES_CGS_REFINE_ALWAYS):
518:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
519:       break;
520:     case (KSP_GMRES_CGS_REFINE_IFNEEDED):
521:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
522:       break;
523:     default:
524:       SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Unknown orthogonalization");
525:     }
526:   } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
527:     cstr = "Modified Gram-Schmidt Orthogonalization";
528:   } else {
529:     cstr = "unknown orthogonalization";
530:   }
531:   if (iascii) {
532:     PetscViewerASCIIPrintf(viewer,"  restart=%D, using %s\n",gmres->max_k,cstr);
533:     PetscViewerASCIIPrintf(viewer,"  happy breakdown tolerance %g\n",(double)gmres->haptol);
534:   } else if (isstring) {
535:     PetscViewerStringSPrintf(viewer,"%s restart %D",cstr,gmres->max_k);
536:   }
537:   return(0);
538: }

540: /*@C
541:    KSPGMRESMonitorKrylov - Calls VecView() for each new direction in the GMRES accumulated Krylov space.

543:    Collective on ksp

545:    Input Parameters:
546: +  ksp - the KSP context
547: .  its - iteration number
548: .  fgnorm - 2-norm of residual (or gradient)
549: -  dummy - an collection of viewers created with KSPViewerCreate()

551:    Options Database Keys:
552: .   -ksp_gmres_kyrlov_monitor

554:    Notes:
555:     A new PETSCVIEWERDRAW is created for each Krylov vector so they can all be simultaneously viewed
556:    Level: intermediate

558: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView(), KSPViewersCreate(), KSPViewersDestroy()
559: @*/
560: PetscErrorCode  KSPGMRESMonitorKrylov(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
561: {
562:   PetscViewers   viewers = (PetscViewers)dummy;
563:   KSP_GMRES      *gmres  = (KSP_GMRES*)ksp->data;
565:   Vec            x;
566:   PetscViewer    viewer;
567:   PetscBool      flg;

570:   PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
571:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&flg);
572:   if (!flg) {
573:     PetscViewerSetType(viewer,PETSCVIEWERDRAW);
574:     PetscViewerDrawSetInfo(viewer,NULL,"Krylov GMRES Monitor",PETSC_DECIDE,PETSC_DECIDE,300,300);
575:   }
576:   x    = VEC_VV(gmres->it+1);
577:   VecView(x,viewer);
578:   return(0);
579: }

581: PetscErrorCode KSPSetFromOptions_GMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
582: {
584:   PetscInt       restart;
585:   PetscReal      haptol,breakdowntol;
586:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
587:   PetscBool      flg;

590:   PetscOptionsHead(PetscOptionsObject,"KSP GMRES Options");
591:   PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
592:   if (flg) { KSPGMRESSetRestart(ksp,restart); }
593:   PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
594:   if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
595:   PetscOptionsReal("-ksp_gmres_breakdown_tolerance","Divergence breakdown tolerance during GMRES restart","KSPGMRESSetBreakdownTolerance",gmres->breakdowntol,&breakdowntol,&flg);
596:   if (flg) { KSPGMRESSetBreakdownTolerance(ksp,breakdowntol); }
597:   flg  = PETSC_FALSE;
598:   PetscOptionsBool("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",flg,&flg,NULL);
599:   if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
600:   PetscOptionsBoolGroupBegin("-ksp_gmres_classicalgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
601:   if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
602:   PetscOptionsBoolGroupEnd("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
603:   if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
604:   PetscOptionsEnum("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType",
605:                           KSPGMRESCGSRefinementTypes,(PetscEnum)gmres->cgstype,(PetscEnum*)&gmres->cgstype,&flg);
606:   flg  = PETSC_FALSE;
607:   PetscOptionsBool("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPMonitorSet",flg,&flg,NULL);
608:   if (flg) {
609:     PetscViewers viewers;
610:     PetscViewersCreate(PetscObjectComm((PetscObject)ksp),&viewers);
611:     KSPMonitorSet(ksp,KSPGMRESMonitorKrylov,viewers,(PetscErrorCode (*)(void**))PetscViewersDestroy);
612:   }
613:   PetscOptionsTail();
614:   return(0);
615: }

617: PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
618: {
619:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

622:   if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
623:   gmres->haptol = tol;
624:   return(0);
625: }

627: PetscErrorCode  KSPGMRESSetBreakdownTolerance_GMRES(KSP ksp,PetscReal tol)
628: {
629:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

632:   if (tol == PETSC_DEFAULT) {
633:     gmres->breakdowntol = 0.1;
634:     return(0);
635:   }
636:   if (tol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Breakdown tolerance must be non-negative");
637:   gmres->breakdowntol = tol;
638:   return(0);
639: }

641: PetscErrorCode  KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
642: {
643:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

646:   *max_k = gmres->max_k;
647:   return(0);
648: }

650: PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
651: {
652:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

656:   if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
657:   if (!ksp->setupstage) {
658:     gmres->max_k = max_k;
659:   } else if (gmres->max_k != max_k) {
660:     gmres->max_k    = max_k;
661:     ksp->setupstage = KSP_SETUP_NEW;
662:     /* free the data structures, then create them again */
663:     KSPReset_GMRES(ksp);
664:   }
665:   return(0);
666: }

668: PetscErrorCode  KSPGMRESSetOrthogonalization_GMRES(KSP ksp,FCN fcn)
669: {
671:   ((KSP_GMRES*)ksp->data)->orthog = fcn;
672:   return(0);
673: }

675: PetscErrorCode  KSPGMRESGetOrthogonalization_GMRES(KSP ksp,FCN *fcn)
676: {
678:   *fcn = ((KSP_GMRES*)ksp->data)->orthog;
679:   return(0);
680: }

682: PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
683: {
684:   KSP_GMRES *gmres;

687:   gmres = (KSP_GMRES*)ksp->data;
688:   gmres->q_preallocate = 1;
689:   return(0);
690: }

692: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
693: {
694:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

697:   gmres->cgstype = type;
698:   return(0);
699: }

701: PetscErrorCode  KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
702: {
703:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

706:   *type = gmres->cgstype;
707:   return(0);
708: }

710: /*@
711:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
712:          in the classical Gram Schmidt orthogonalization.

714:    Logically Collective on ksp

716:    Input Parameters:
717: +  ksp - the Krylov space context
718: -  type - the type of refinement

720:   Options Database:
721: .  -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>

723:    Level: intermediate

725: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
726:           KSPGMRESGetOrthogonalization()
727: @*/
728: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
729: {

735:   PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
736:   return(0);
737: }

739: /*@
740:    KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
741:          in the classical Gram Schmidt orthogonalization.

743:    Not Collective

745:    Input Parameter:
746: .  ksp - the Krylov space context

748:    Output Parameter:
749: .  type - the type of refinement

751:   Options Database:
752: .  -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>

754:    Level: intermediate

756: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
757:           KSPGMRESGetOrthogonalization()
758: @*/
759: PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
760: {

765:   PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
766:   return(0);
767: }


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

773:    Logically Collective on ksp

775:    Input Parameters:
776: +  ksp - the Krylov space context
777: -  restart - integer restart value

779:   Options Database:
780: .  -ksp_gmres_restart <positive integer>

782:     Note: The default value is 30.

784:    Level: intermediate

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: }

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

802:    Not Collective

804:    Input Parameter:
805: .  ksp - the Krylov space context

807:    Output Parameter:
808: .   restart - integer restart value

810:     Note: The default value is 30.

812:    Level: intermediate

814: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetRestart()
815: @*/
816: PetscErrorCode  KSPGMRESGetRestart(KSP ksp, PetscInt *restart)
817: {

821:   PetscUseMethod(ksp,"KSPGMRESGetRestart_C",(KSP,PetscInt*),(ksp,restart));
822:   return(0);
823: }

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

828:    Logically Collective on ksp

830:    Input Parameters:
831: +  ksp - the Krylov space context
832: -  tol - the tolerance

834:   Options Database:
835: .  -ksp_gmres_haptol <positive real value>

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

841:    Level: intermediate

843: .seealso: KSPSetTolerances()
844: @*/
845: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
846: {

851:   PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
852:   return(0);
853: }

855: /*@
856:    KSPGMRESSetBreakdownTolerance - Sets tolerance for determining divergence breakdown in GMRES.

858:    Logically Collective on ksp

860:    Input Parameters:
861: +  ksp - the Krylov space context
862: -  tol - the tolerance

864:   Options Database:
865: .  -ksp_gmres_breakdown_tolerance <positive real value>

867:    Note: divergence breakdown occurs when GMRES residual increases significantly
868:          during restart

870:    Level: intermediate

872: .seealso: KSPSetTolerances(), KSPGMRESSetHapTol()
873: @*/
874: PetscErrorCode  KSPGMRESSetBreakdownTolerance(KSP ksp,PetscReal tol)
875: {

880:   PetscTryMethod((ksp),"KSPGMRESSetBreakdownTolerance_C",(KSP,PetscReal),(ksp,tol));
881:   return(0);
882: }

884: /*MC
885:      KSPGMRES - Implements the Generalized Minimal Residual method.
886:                 (Saad and Schultz, 1986) with restart


889:    Options Database Keys:
890: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
891: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
892: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
893:                              vectors are allocated as needed)
894: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
895: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
896: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
897:                                    stability of the classical Gram-Schmidt  orthogonalization.
898: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

900:    Level: beginner

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

905:    References:
906: .     1. - YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
907:           SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.

909: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPLGMRES,
910:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
911:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
912:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()

914: M*/

916: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
917: {
918:   KSP_GMRES      *gmres;

922:   PetscNewLog(ksp,&gmres);
923:   ksp->data = (void*)gmres;

925:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,4);
926:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
927:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,2);
928:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
929:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

931:   ksp->ops->buildsolution                = KSPBuildSolution_GMRES;
932:   ksp->ops->setup                        = KSPSetUp_GMRES;
933:   ksp->ops->solve                        = KSPSolve_GMRES;
934:   ksp->ops->reset                        = KSPReset_GMRES;
935:   ksp->ops->destroy                      = KSPDestroy_GMRES;
936:   ksp->ops->view                         = KSPView_GMRES;
937:   ksp->ops->setfromoptions               = KSPSetFromOptions_GMRES;
938:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
939:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;
940: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_ESSL)
941:   ksp->ops->computeritz                  = KSPComputeRitz_GMRES;
942: #endif
943:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
944:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
945:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
946:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
947:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
948:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
949:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetBreakdownTolerance_C",KSPGMRESSetBreakdownTolerance_GMRES);
950:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
951:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

953:   gmres->haptol         = 1.0e-30;
954:   gmres->breakdowntol   = 0.1;
955:   gmres->q_preallocate  = 0;
956:   gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
957:   gmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
958:   gmres->nrs            = NULL;
959:   gmres->sol_temp       = NULL;
960:   gmres->max_k          = GMRES_DEFAULT_MAXK;
961:   gmres->Rsvd           = NULL;
962:   gmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
963:   gmres->orthogwork     = NULL;
964:   return(0);
965: }