Actual source code: gmres.c

petsc-3.15.0 2021-04-05
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:   KSPLogErrorHistory(ksp);
145:   KSPMonitor(ksp,ksp->its,res);
146:   if (!res) {
147:     ksp->reason = KSP_CONVERGED_ATOL;
148:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
149:     return(0);
150:   }

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

165:     /* update hessenberg matrix and do Gram-Schmidt */
166:     (*gmres->orthog)(ksp,it);
167:     if (ksp->reason) break;

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

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

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

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

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

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

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

215:   if (itcount) *itcount = it;


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

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

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

239:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
240:   ksp->its = 0;
241:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

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

277: PetscErrorCode KSPReset_GMRES(KSP ksp)
278: {
279:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
281:   PetscInt       i;

284:   /* Free the Hessenberg matrices */
285:   PetscFree5(gmres->hh_origin,gmres->hes_origin,gmres->rs_origin,gmres->cc_origin,gmres->ss_origin);
286:   PetscFree(gmres->hes_ritz);

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

298:   PetscFree(gmres->user_work);
299:   PetscFree(gmres->mwork_alloc);
300:   PetscFree(gmres->nrs);
301:   VecDestroy(&gmres->sol_temp);
302:   PetscFree(gmres->Rsvd);
303:   PetscFree(gmres->Dsvd);
304:   PetscFree(gmres->orthogwork);

306:   gmres->vv_allocated   = 0;
307:   gmres->vecs_allocated = 0;
308:   gmres->sol_temp       = NULL;
309:   return(0);
310: }

312: PetscErrorCode KSPDestroy_GMRES(KSP ksp)
313: {

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

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

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

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

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

364:     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)));
365:     return(0);
366:   }
367:   for (ii=1; ii<=it; ii++) {
368:     k  = it - ii;
369:     tt = *GRS(k);
370:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
371:     if (*HH(k,k) == 0.0) {
372:       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);
373:       else {
374:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
375:         PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; k = %D\n",k);
376:         return(0);
377:       }
378:     }
379:     nrs[k] = tt / *HH(k,k);
380:   }

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

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

404:   hh = HH(0,it);
405:   cc = CC(0);
406:   ss = SS(0);

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

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

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

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

468:   gmres->vv_allocated += nalloc;

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

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

481: PetscErrorCode KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
482: {
483:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

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

500:   KSPGMRESBuildSoln(gmres->nrs,ksp->vec_sol,ptr,ksp,gmres->it);
501:   if (result) *result = ptr;
502:   return(0);
503: }

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

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

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

546:    Collective on ksp

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

554:    Options Database Keys:
555: .   -ksp_gmres_kyrlov_monitor

557:    Notes:
558:     A new PETSCVIEWERDRAW is created for each Krylov vector so they can all be simultaneously viewed
559:    Level: intermediate

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

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

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

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

620: PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP ksp,PetscReal tol)
621: {
622:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

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

630: PetscErrorCode  KSPGMRESSetBreakdownTolerance_GMRES(KSP ksp,PetscReal tol)
631: {
632:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

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

644: PetscErrorCode  KSPGMRESGetRestart_GMRES(KSP ksp,PetscInt *max_k)
645: {
646:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

649:   *max_k = gmres->max_k;
650:   return(0);
651: }

653: PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP ksp,PetscInt max_k)
654: {
655:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;

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

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

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

685: PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
686: {
687:   KSP_GMRES *gmres;

690:   gmres = (KSP_GMRES*)ksp->data;
691:   gmres->q_preallocate = 1;
692:   return(0);
693: }

695: PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType type)
696: {
697:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

700:   gmres->cgstype = type;
701:   return(0);
702: }

704: PetscErrorCode  KSPGMRESGetCGSRefinementType_GMRES(KSP ksp,KSPGMRESCGSRefinementType *type)
705: {
706:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;

709:   *type = gmres->cgstype;
710:   return(0);
711: }

713: /*@
714:    KSPGMRESSetCGSRefinementType - Sets the type of iterative refinement to use
715:          in the classical Gram Schmidt orthogonalization.

717:    Logically Collective on ksp

719:    Input Parameters:
720: +  ksp - the Krylov space context
721: -  type - the type of refinement

723:   Options Database:
724: .  -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>

726:    Level: intermediate

728: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType(),
729:           KSPGMRESGetOrthogonalization()
730: @*/
731: PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
732: {

738:   PetscTryMethod(ksp,"KSPGMRESSetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType),(ksp,type));
739:   return(0);
740: }

742: /*@
743:    KSPGMRESGetCGSRefinementType - Gets the type of iterative refinement to use
744:          in the classical Gram Schmidt orthogonalization.

746:    Not Collective

748:    Input Parameter:
749: .  ksp - the Krylov space context

751:    Output Parameter:
752: .  type - the type of refinement

754:   Options Database:
755: .  -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>

757:    Level: intermediate

759: .seealso: KSPGMRESSetOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetCGSRefinementType(),
760:           KSPGMRESGetOrthogonalization()
761: @*/
762: PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType *type)
763: {

768:   PetscUseMethod(ksp,"KSPGMRESGetCGSRefinementType_C",(KSP,KSPGMRESCGSRefinementType*),(ksp,type));
769:   return(0);
770: }


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

776:    Logically Collective on ksp

778:    Input Parameters:
779: +  ksp - the Krylov space context
780: -  restart - integer restart value

782:   Options Database:
783: .  -ksp_gmres_restart <positive integer>

785:     Note: The default value is 30.

787:    Level: intermediate

789: .seealso: KSPSetTolerances(), KSPGMRESSetOrthogonalization(), KSPGMRESSetPreAllocateVectors(), KSPGMRESGetRestart()
790: @*/
791: PetscErrorCode  KSPGMRESSetRestart(KSP ksp, PetscInt restart)
792: {


798:   PetscTryMethod(ksp,"KSPGMRESSetRestart_C",(KSP,PetscInt),(ksp,restart));
799:   return(0);
800: }

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

805:    Not Collective

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

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

813:     Note: The default value is 30.

815:    Level: intermediate

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

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

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

831:    Logically Collective on ksp

833:    Input Parameters:
834: +  ksp - the Krylov space context
835: -  tol - the tolerance

837:   Options Database:
838: .  -ksp_gmres_haptol <positive real value>

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

844:    Level: intermediate

846: .seealso: KSPSetTolerances()
847: @*/
848: PetscErrorCode  KSPGMRESSetHapTol(KSP ksp,PetscReal tol)
849: {

854:   PetscTryMethod((ksp),"KSPGMRESSetHapTol_C",(KSP,PetscReal),((ksp),(tol)));
855:   return(0);
856: }

858: /*@
859:    KSPGMRESSetBreakdownTolerance - Sets tolerance for determining divergence breakdown in GMRES.

861:    Logically Collective on ksp

863:    Input Parameters:
864: +  ksp - the Krylov space context
865: -  tol - the tolerance

867:   Options Database:
868: .  -ksp_gmres_breakdown_tolerance <positive real value>

870:    Note: divergence breakdown occurs when GMRES residual increases significantly
871:          during restart

873:    Level: intermediate

875: .seealso: KSPSetTolerances(), KSPGMRESSetHapTol()
876: @*/
877: PetscErrorCode  KSPGMRESSetBreakdownTolerance(KSP ksp,PetscReal tol)
878: {

883:   PetscTryMethod((ksp),"KSPGMRESSetBreakdownTolerance_C",(KSP,PetscReal),(ksp,tol));
884:   return(0);
885: }

887: /*MC
888:      KSPGMRES - Implements the Generalized Minimal Residual method.
889:                 (Saad and Schultz, 1986) with restart


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

903:    Level: beginner

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

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

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

917: M*/

919: PETSC_EXTERN PetscErrorCode KSPCreate_GMRES(KSP ksp)
920: {
921:   KSP_GMRES      *gmres;

925:   PetscNewLog(ksp,&gmres);
926:   ksp->data = (void*)gmres;

928:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,4);
929:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
930:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,2);
931:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
932:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

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

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