Actual source code: fgmres.c


  2: /*
  3:     This file implements FGMRES (a Generalized Minimal Residual) method.
  4:     Reference:  Saad, 1993.

  6:     Preconditioning:  If the preconditioner is constant then this fgmres
  7:     code is equivalent to RIGHT-PRECONDITIONED GMRES.
  8:     FGMRES is a modification of gmres that allows the preconditioner to change
  9:     at each iteration.

 11:     Restarts:  Restarts are basically solves with x0 not equal to zero.

 13:        Contributed by Allison Baker

 15: */

 17: #include <../src/ksp/ksp/impls/gmres/fgmres/fgmresimpl.h>
 18: #define FGMRES_DELTA_DIRECTIONS 10
 19: #define FGMRES_DEFAULT_MAXK     30
 20: static PetscErrorCode KSPFGMRESGetNewVectors(KSP, PetscInt);
 21: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP, PetscInt, PetscBool, PetscReal *);
 22: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);

 24: /*

 26:     KSPSetUp_FGMRES - Sets up the workspace needed by fgmres.

 28:     This is called once, usually automatically by KSPSolve() or KSPSetUp(),
 29:     but can be called directly by KSPSetUp().

 31: */
 32: PetscErrorCode KSPSetUp_FGMRES(KSP ksp)
 33: {
 34:   PetscInt    max_k, k;
 35:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;

 37:   max_k = fgmres->max_k;

 39:   KSPSetUp_GMRES(ksp);

 41:   PetscMalloc1(max_k + 2, &fgmres->prevecs);
 42:   PetscMalloc1(max_k + 2, &fgmres->prevecs_user_work);

 44:   /* fgmres->vv_allocated includes extra work vectors, which are not used in the additional
 45:      block of vectors used to store the preconditioned directions, hence  the -VEC_OFFSET
 46:      term for this first allocation of vectors holding preconditioned directions */
 47:   KSPCreateVecs(ksp, fgmres->vv_allocated - VEC_OFFSET, &fgmres->prevecs_user_work[0], 0, NULL);
 48:   for (k = 0; k < fgmres->vv_allocated - VEC_OFFSET; k++) fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
 49:   return 0;
 50: }

 52: /*
 53:     KSPFGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED)
 54: */
 55: static PetscErrorCode KSPFGMRESResidual(KSP ksp)
 56: {
 57:   KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
 58:   Mat         Amat, Pmat;

 60:   PCGetOperators(ksp->pc, &Amat, &Pmat);

 62:   /* put A*x into VEC_TEMP */
 63:   KSP_MatMult(ksp, Amat, ksp->vec_sol, VEC_TEMP);
 64:   /* now put residual (-A*x + f) into vec_vv(0) */
 65:   VecWAXPY(VEC_VV(0), -1.0, VEC_TEMP, ksp->vec_rhs);
 66:   return 0;
 67: }

 69: /*

 71:     KSPFGMRESCycle - Run fgmres, possibly with restart.  Return residual
 72:                   history if requested.

 74:     input parameters:
 75: .        fgmres  - structure containing parameters and work areas

 77:     output parameters:
 78: .        itcount - number of iterations used.  If null, ignored.
 79: .        converged - 0 if not converged

 81:     Notes:
 82:     On entry, the value in vector VEC_VV(0) should be
 83:     the initial residual.

 85:  */
 86: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount, KSP ksp)
 87: {
 88:   KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
 89:   PetscReal   res_norm;
 90:   PetscReal   hapbnd, tt;
 91:   PetscBool   hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
 92:   PetscInt    loc_it;                /* local count of # of dir. in Krylov space */
 93:   PetscInt    max_k = fgmres->max_k; /* max # of directions Krylov space */
 94:   Mat         Amat, Pmat;

 96:   /* Number of pseudo iterations since last restart is the number
 97:      of prestart directions */
 98:   loc_it = 0;

100:   /* note: (fgmres->it) is always set one less than (loc_it) It is used in
101:      KSPBUILDSolution_FGMRES, where it is passed to KSPFGMRESBuildSoln.
102:      Note that when KSPFGMRESBuildSoln is called from this function,
103:      (loc_it -1) is passed, so the two are equivalent */
104:   fgmres->it = (loc_it - 1);

106:   /* initial residual is in VEC_VV(0)  - compute its norm*/
107:   VecNorm(VEC_VV(0), NORM_2, &res_norm);
108:   KSPCheckNorm(ksp, res_norm);

110:   /* first entry in right-hand-side of hessenberg system is just
111:      the initial residual norm */
112:   *RS(0) = res_norm;

114:   ksp->rnorm = res_norm;
115:   KSPLogResidualHistory(ksp, res_norm);
116:   KSPMonitor(ksp, ksp->its, res_norm);

118:   /* check for the convergence - maybe the current guess is good enough */
119:   (*ksp->converged)(ksp, ksp->its, res_norm, &ksp->reason, ksp->cnvP);
120:   if (ksp->reason) {
121:     if (itcount) *itcount = 0;
122:     return 0;
123:   }

125:   /* scale VEC_VV (the initial residual) */
126:   VecScale(VEC_VV(0), 1.0 / res_norm);

128:   /* MAIN ITERATION LOOP BEGINNING*/
129:   /* keep iterating until we have converged OR generated the max number
130:      of directions OR reached the max number of iterations for the method */
131:   while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
132:     if (loc_it) {
133:       KSPLogResidualHistory(ksp, res_norm);
134:       KSPMonitor(ksp, ksp->its, res_norm);
135:     }
136:     fgmres->it = (loc_it - 1);

138:     /* see if more space is needed for work vectors */
139:     if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
140:       KSPFGMRESGetNewVectors(ksp, loc_it + 1);
141:       /* (loc_it+1) is passed in as number of the first vector that should
142:          be allocated */
143:     }

145:     /* CHANGE THE PRECONDITIONER? */
146:     /* ModifyPC is the callback function that can be used to
147:        change the PC or its attributes before its applied */
148:     (*fgmres->modifypc)(ksp, ksp->its, loc_it, res_norm, fgmres->modifyctx);

150:     /* apply PRECONDITIONER to direction vector and store with
151:        preconditioned vectors in prevec */
152:     KSP_PCApply(ksp, VEC_VV(loc_it), PREVEC(loc_it));

154:     PCGetOperators(ksp->pc, &Amat, &Pmat);
155:     /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
156:     KSP_MatMult(ksp, Amat, PREVEC(loc_it), VEC_VV(1 + loc_it));

158:     /* update hessenberg matrix and do Gram-Schmidt - new direction is in
159:        VEC_VV(1+loc_it)*/
160:     (*fgmres->orthog)(ksp, loc_it);

162:     /* new entry in hessenburg is the 2-norm of our new direction */
163:     VecNorm(VEC_VV(loc_it + 1), NORM_2, &tt);
164:     KSPCheckNorm(ksp, tt);

166:     *HH(loc_it + 1, loc_it)  = tt;
167:     *HES(loc_it + 1, loc_it) = tt;

169:     /* Happy Breakdown Check */
170:     hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
171:     /* RS(loc_it) contains the res_norm from the last iteration  */
172:     hapbnd = PetscMin(fgmres->haptol, hapbnd);
173:     if (tt > hapbnd) {
174:       /* scale new direction by its norm */
175:       VecScale(VEC_VV(loc_it + 1), 1.0 / tt);
176:     } else {
177:       /* This happens when the solution is exactly reached. */
178:       /* So there is no new direction... */
179:       VecSet(VEC_TEMP, 0.0); /* set VEC_TEMP to 0 */
180:       hapend = PETSC_TRUE;
181:     }
182:     /* note that for FGMRES we could get HES(loc_it+1, loc_it)  = 0 and the
183:        current solution would not be exact if HES was singular.  Note that
184:        HH non-singular implies that HES is no singular, and HES is guaranteed
185:        to be nonsingular when PREVECS are linearly independent and A is
186:        nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
187:        of HES). So we should really add a check to verify that HES is nonsingular.*/

189:     /* Now apply rotations to new col of hessenberg (and right side of system),
190:        calculate new rotation, and get new residual norm at the same time*/
191:     KSPFGMRESUpdateHessenberg(ksp, loc_it, hapend, &res_norm);
192:     if (ksp->reason) break;

194:     loc_it++;
195:     fgmres->it = (loc_it - 1); /* Add this here in case it has converged */

197:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
198:     ksp->its++;
199:     ksp->rnorm = res_norm;
200:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

202:     (*ksp->converged)(ksp, ksp->its, res_norm, &ksp->reason, ksp->cnvP);

204:     /* Catch error in happy breakdown and signal convergence and break from loop */
205:     if (hapend) {
206:       if (!ksp->reason) {
208:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
209:         break;
210:       }
211:     }
212:   }
213:   /* END OF ITERATION LOOP */
214:   KSPLogResidualHistory(ksp, res_norm);

216:   /*
217:      Monitor if we know that we will not return for a restart */
218:   if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) KSPMonitor(ksp, ksp->its, res_norm);

220:   if (itcount) *itcount = loc_it;

222:   /*
223:     Down here we have to solve for the "best" coefficients of the Krylov
224:     columns, add the solution values together, and possibly unwind the
225:     preconditioning from the solution
226:    */

228:   /* Form the solution (or the solution so far) */
229:   /* Note: must pass in (loc_it-1) for iteration count so that KSPFGMRESBuildSoln
230:      properly navigates */

232:   KSPFGMRESBuildSoln(RS(0), ksp->vec_sol, ksp->vec_sol, ksp, loc_it - 1);
233:   return 0;
234: }

236: /*
237:     KSPSolve_FGMRES - This routine applies the FGMRES method.

239:    Input Parameter:
240: .     ksp - the Krylov space object that was set to use fgmres

242:    Output Parameter:
243: .     outits - number of iterations used

245: */

247: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
248: {
249:   PetscInt    cycle_its = 0; /* iterations done in a call to KSPFGMRESCycle */
250:   KSP_FGMRES *fgmres    = (KSP_FGMRES *)ksp->data;
251:   PetscBool   diagonalscale;

253:   PCGetDiagonalScale(ksp->pc, &diagonalscale);

256:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
257:   ksp->its = 0;
258:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

260:   /* Compute the initial (NOT preconditioned) residual */
261:   if (!ksp->guess_zero) {
262:     KSPFGMRESResidual(ksp);
263:   } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
264:     VecCopy(ksp->vec_rhs, VEC_VV(0));
265:   }
266:   /* This may be true only on a subset of MPI ranks; setting it here so it will be detected by the first norm computaion in the Krylov method */
267:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) VecSetInf(VEC_VV(0));

269:   /* now the residual is in VEC_VV(0) - which is what
270:      KSPFGMRESCycle expects... */

272:   KSPFGMRESCycle(&cycle_its, ksp);
273:   while (!ksp->reason) {
274:     KSPFGMRESResidual(ksp);
275:     if (ksp->its >= ksp->max_it) break;
276:     KSPFGMRESCycle(&cycle_its, ksp);
277:   }
278:   /* mark lack of convergence */
279:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
280:   return 0;
281: }

283: extern PetscErrorCode KSPReset_FGMRES(KSP);
284: /*

286:    KSPDestroy_FGMRES - Frees all memory space used by the Krylov method.

288: */
289: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
290: {
291:   KSPReset_FGMRES(ksp);
292:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFGMRESSetModifyPC_C", NULL);
293:   KSPDestroy_GMRES(ksp);
294:   return 0;
295: }

297: /*
298:     KSPFGMRESBuildSoln - create the solution from the starting vector and the
299:                       current iterates.

301:     Input parameters:
302:         nrs - work area of size it + 1.
303:         vguess  - index of initial guess
304:         vdest - index of result.  Note that vguess may == vdest (replace
305:                 guess with the solution).
306:         it - HH upper triangular part is a block of size (it+1) x (it+1)

308:      This is an internal routine that knows about the FGMRES internals.
309:  */
310: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar *nrs, Vec vguess, Vec vdest, KSP ksp, PetscInt it)
311: {
312:   PetscScalar tt;
313:   PetscInt    ii, k, j;
314:   KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);

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

318:   /* If it is < 0, no fgmres steps have been performed */
319:   if (it < 0) {
320:     VecCopy(vguess, vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
321:     return 0;
322:   }

324:   /* so fgmres steps HAVE been performed */

326:   /* solve the upper triangular system - RS is the right side and HH is
327:      the upper triangular matrix  - put soln in nrs */
328:   if (*HH(it, it) != 0.0) {
329:     nrs[it] = *RS(it) / *HH(it, it);
330:   } else {
331:     nrs[it] = 0.0;
332:   }
333:   for (ii = 1; ii <= it; ii++) {
334:     k  = it - ii;
335:     tt = *RS(k);
336:     for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
337:     nrs[k] = tt / *HH(k, k);
338:   }

340:   /* Accumulate the correction to the soln of the preconditioned prob. in
341:      VEC_TEMP - note that we use the preconditioned vectors  */
342:   VecSet(VEC_TEMP, 0.0); /* set VEC_TEMP components to 0 */
343:   VecMAXPY(VEC_TEMP, it + 1, nrs, &PREVEC(0));

345:   /* put updated solution into vdest.*/
346:   if (vdest != vguess) {
347:     VecCopy(VEC_TEMP, vdest);
348:     VecAXPY(vdest, 1.0, vguess);
349:   } else { /* replace guess with solution */
350:     VecAXPY(vdest, 1.0, VEC_TEMP);
351:   }
352:   return 0;
353: }

355: /*

357:     KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
358:                             Return new residual.

360:     input parameters:

362: .        ksp -    Krylov space object
363: .        it  -    plane rotations are applied to the (it+1)th column of the
364:                   modified hessenberg (i.e. HH(:,it))
365: .        hapend - PETSC_FALSE not happy breakdown ending.

367:     output parameters:
368: .        res - the new residual

370:  */
371: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
372: {
373:   PetscScalar *hh, *cc, *ss, tt;
374:   PetscInt     j;
375:   KSP_FGMRES  *fgmres = (KSP_FGMRES *)(ksp->data);

377:   hh = HH(0, it); /* pointer to beginning of column to update - so
378:                       incrementing hh "steps down" the (it+1)th col of HH*/
379:   cc = CC(0);     /* beginning of cosine rotations */
380:   ss = SS(0);     /* beginning of sine rotations */

382:   /* Apply all the previously computed plane rotations to the new column
383:      of the Hessenberg matrix */
384:   /* Note: this uses the rotation [conj(c)  s ; -s   c], c= cos(theta), s= sin(theta),
385:      and some refs have [c   s ; -conj(s)  c] (don't be confused!) */

387:   for (j = 1; j <= it; j++) {
388:     tt  = *hh;
389:     *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
390:     hh++;
391:     *hh = *cc++ * *hh - (*ss++ * tt);
392:     /* hh, cc, and ss have all been incremented one by end of loop */
393:   }

395:   /*
396:     compute the new plane rotation, and apply it to:
397:      1) the right-hand-side of the Hessenberg system (RS)
398:         note: it affects RS(it) and RS(it+1)
399:      2) the new column of the Hessenberg matrix
400:         note: it affects HH(it,it) which is currently pointed to
401:         by hh and HH(it+1, it) (*(hh+1))
402:     thus obtaining the updated value of the residual...
403:   */

405:   /* compute new plane rotation */

407:   if (!hapend) {
408:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
409:     if (tt == 0.0) {
410:       ksp->reason = KSP_DIVERGED_NULL;
411:       return 0;
412:     }

414:     *cc = *hh / tt;       /* new cosine value */
415:     *ss = *(hh + 1) / tt; /* new sine value */

417:     /* apply to 1) and 2) */
418:     *RS(it + 1) = -(*ss * *RS(it));
419:     *RS(it)     = PetscConj(*cc) * *RS(it);
420:     *hh         = PetscConj(*cc) * *hh + *ss * *(hh + 1);

422:     /* residual is the last element (it+1) of right-hand side! */
423:     *res = PetscAbsScalar(*RS(it + 1));

425:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
426:             another rotation matrix (so RH doesn't change).  The new residual is
427:             always the new sine term times the residual from last time (RS(it)),
428:             but now the new sine rotation would be zero...so the residual should
429:             be zero...so we will multiply "zero" by the last residual.  This might
430:             not be exactly what we want to do here -could just return "zero". */

432:     *res = 0.0;
433:   }
434:   return 0;
435: }

437: /*

439:    KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
440:                          VEC_VV(it), and more preconditioned work vectors, starting
441:                          from PREVEC(i).

443: */
444: static PetscErrorCode KSPFGMRESGetNewVectors(KSP ksp, PetscInt it)
445: {
446:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
447:   PetscInt    nwork  = fgmres->nwork_alloc; /* number of work vector chunks allocated */
448:   PetscInt    nalloc;                       /* number to allocate */
449:   PetscInt    k;

451:   nalloc = fgmres->delta_allocate; /* number of vectors to allocate
452:                                       in a single chunk */

454:   /* Adjust the number to allocate to make sure that we don't exceed the
455:      number of available slots (fgmres->vecs_allocated)*/
456:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated) nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
457:   if (!nalloc) return 0;

459:   fgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */

461:   /* work vectors */
462:   KSPCreateVecs(ksp, nalloc, &fgmres->user_work[nwork], 0, NULL);
463:   for (k = 0; k < nalloc; k++) fgmres->vecs[it + VEC_OFFSET + k] = fgmres->user_work[nwork][k];
464:   /* specify size of chunk allocated */
465:   fgmres->mwork_alloc[nwork] = nalloc;

467:   /* preconditioned vectors */
468:   KSPCreateVecs(ksp, nalloc, &fgmres->prevecs_user_work[nwork], 0, NULL);
469:   for (k = 0; k < nalloc; k++) fgmres->prevecs[it + k] = fgmres->prevecs_user_work[nwork][k];

471:   /* increment the number of work vector chunks */
472:   fgmres->nwork_alloc++;
473:   return 0;
474: }

476: /*

478:    KSPBuildSolution_FGMRES

480:      Input Parameter:
481: .     ksp - the Krylov space object
482: .     ptr-

484:    Output Parameter:
485: .     result - the solution

487:    Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
488:    calls directly.

490: */
491: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp, Vec ptr, Vec *result)
492: {
493:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;

495:   if (!ptr) {
496:     if (!fgmres->sol_temp) { VecDuplicate(ksp->vec_sol, &fgmres->sol_temp); }
497:     ptr = fgmres->sol_temp;
498:   }
499:   if (!fgmres->nrs) {
500:     /* allocate the work area */
501:     PetscMalloc1(fgmres->max_k, &fgmres->nrs);
502:   }

504:   KSPFGMRESBuildSoln(fgmres->nrs, ksp->vec_sol, ptr, ksp, fgmres->it);
505:   if (result) *result = ptr;
506:   return 0;
507: }

509: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
510: {
511:   PetscBool flg;

513:   KSPSetFromOptions_GMRES(ksp, PetscOptionsObject);
514:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP flexible GMRES Options");
515:   PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange", "do not vary the preconditioner", "KSPFGMRESSetModifyPC", &flg);
516:   if (flg) KSPFGMRESSetModifyPC(ksp, KSPFGMRESModifyPCNoChange, NULL, NULL);
517:   PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp", "vary the KSP based preconditioner", "KSPFGMRESSetModifyPC", &flg);
518:   if (flg) KSPFGMRESSetModifyPC(ksp, KSPFGMRESModifyPCKSP, NULL, NULL);
519:   PetscOptionsHeadEnd();
520:   return 0;
521: }

523: typedef PetscErrorCode (*FCN1)(KSP, PetscInt, PetscInt, PetscReal, void *); /* force argument to next function to not be extern C*/
524: typedef PetscErrorCode (*FCN2)(void *);

526: static PetscErrorCode KSPFGMRESSetModifyPC_FGMRES(KSP ksp, FCN1 fcn, void *ctx, FCN2 d)
527: {
529:   ((KSP_FGMRES *)ksp->data)->modifypc      = fcn;
530:   ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
531:   ((KSP_FGMRES *)ksp->data)->modifyctx     = ctx;
532:   return 0;
533: }

535: PetscErrorCode KSPReset_FGMRES(KSP ksp)
536: {
537:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
538:   PetscInt    i;

540:   PetscFree(fgmres->prevecs);
541:   if (fgmres->nwork_alloc > 0) {
542:     i = 0;
543:     /* In the first allocation we allocated VEC_OFFSET fewer vectors in prevecs */
544:     VecDestroyVecs(fgmres->mwork_alloc[i] - VEC_OFFSET, &fgmres->prevecs_user_work[i]);
545:     for (i = 1; i < fgmres->nwork_alloc; i++) VecDestroyVecs(fgmres->mwork_alloc[i], &fgmres->prevecs_user_work[i]);
546:   }
547:   PetscFree(fgmres->prevecs_user_work);
548:   if (fgmres->modifydestroy) (*fgmres->modifydestroy)(fgmres->modifyctx);
549:   KSPReset_GMRES(ksp);
550:   return 0;
551: }

553: PetscErrorCode KSPGMRESSetRestart_FGMRES(KSP ksp, PetscInt max_k)
554: {
555:   KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;

558:   if (!ksp->setupstage) {
559:     gmres->max_k = max_k;
560:   } else if (gmres->max_k != max_k) {
561:     gmres->max_k    = max_k;
562:     ksp->setupstage = KSP_SETUP_NEW;
563:     /* free the data structures, then create them again */
564:     KSPReset_FGMRES(ksp);
565:   }
566:   return 0;
567: }

569: PetscErrorCode KSPGMRESGetRestart_FGMRES(KSP ksp, PetscInt *max_k)
570: {
571:   KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;

573:   *max_k = gmres->max_k;
574:   return 0;
575: }

577: /*MC
578:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
579:                 developed by Saad with restart

581:    Options Database Keys:
582: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
583: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
584: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
585:                              vectors are allocated as needed)
586: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
587: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
588: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
589:                                    stability of the classical Gram-Schmidt  orthogonalization.
590: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
591: .   -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
592: -   -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

594:    Level: beginner

596:     Notes:
597:     See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
598:            Only right preconditioning is supported.

600:     Notes:
601:     The following options -ksp_type fgmres -pc_type ksp -ksp_ksp_type bcgs -ksp_view -ksp_pc_type jacobi make the preconditioner (or inner solver)
602:            be bi-CG-stab with a preconditioner of Jacobi.

604:     Developer Notes:
605:     This object is subclassed off of KSPGMRES

607: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPGMRES`, `KSPLGMRES`,
608:           `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
609:           `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
610:           `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPFGMRESSetModifyPC()`,
611:           `KSPFGMRESModifyPCKSP()`

613: M*/

615: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
616: {
617:   KSP_FGMRES *fgmres;

619:   PetscNew(&fgmres);

621:   ksp->data                              = (void *)fgmres;
622:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
623:   ksp->ops->setup                        = KSPSetUp_FGMRES;
624:   ksp->ops->solve                        = KSPSolve_FGMRES;
625:   ksp->ops->reset                        = KSPReset_FGMRES;
626:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
627:   ksp->ops->view                         = KSPView_GMRES;
628:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
629:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
630:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

632:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 3);
633:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

635:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES);
636:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES);
637:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", KSPGMRESGetOrthogonalization_GMRES);
638:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_FGMRES);
639:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_FGMRES);
640:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFGMRESSetModifyPC_C", KSPFGMRESSetModifyPC_FGMRES);
641:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES);
642:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", KSPGMRESGetCGSRefinementType_GMRES);

644:   fgmres->haptol         = 1.0e-30;
645:   fgmres->q_preallocate  = 0;
646:   fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
647:   fgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
648:   fgmres->nrs            = NULL;
649:   fgmres->sol_temp       = NULL;
650:   fgmres->max_k          = FGMRES_DEFAULT_MAXK;
651:   fgmres->Rsvd           = NULL;
652:   fgmres->orthogwork     = NULL;
653:   fgmres->modifypc       = KSPFGMRESModifyPCNoChange;
654:   fgmres->modifyctx      = NULL;
655:   fgmres->modifydestroy  = NULL;
656:   fgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
657:   return 0;
658: }