Actual source code: fgmres.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  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: {
 35:   PetscInt       max_k,k;
 36:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

 39:   max_k = fgmres->max_k;

 41:   KSPSetUp_GMRES(ksp);

 43:   PetscMalloc1(max_k+2,&fgmres->prevecs);
 44:   PetscMalloc1(max_k+2,&fgmres->prevecs_user_work);
 45:   PetscLogObjectMemory((PetscObject)ksp,(max_k+2)*(2*sizeof(void*)));

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

 58: /*
 59:     KSPFGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED)
 60: */
 61: static PetscErrorCode KSPFGMRESResidual(KSP ksp)
 62: {
 63:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)(ksp->data);
 64:   Mat            Amat,Pmat;

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

 70:   /* put A*x into VEC_TEMP */
 71:   KSP_MatMult(ksp,Amat,ksp->vec_sol,VEC_TEMP);
 72:   /* now put residual (-A*x + f) into vec_vv(0) */
 73:   VecWAXPY(VEC_VV(0),-1.0,VEC_TEMP,ksp->vec_rhs);
 74:   return(0);
 75: }

 77: /*

 79:     KSPFGMRESCycle - Run fgmres, possibly with restart.  Return residual
 80:                   history if requested.

 82:     input parameters:
 83: .        fgmres  - structure containing parameters and work areas

 85:     output parameters:
 86: .        itcount - number of iterations used.  If null, ignored.
 87: .        converged - 0 if not converged


 90:     Notes:
 91:     On entry, the value in vector VEC_VV(0) should be
 92:     the initial residual.


 95:  */
 96: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount,KSP ksp)
 97: {

 99:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)(ksp->data);
100:   PetscReal      res_norm;
101:   PetscReal      hapbnd,tt;
102:   PetscBool      hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
104:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
105:   PetscInt       max_k = fgmres->max_k; /* max # of directions Krylov space */
106:   Mat            Amat,Pmat;

109:   /* Number of pseudo iterations since last restart is the number
110:      of prestart directions */
111:   loc_it = 0;

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

119:   /* initial residual is in VEC_VV(0)  - compute its norm*/
120:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
121:   KSPCheckNorm(ksp,res_norm);

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

127:   ksp->rnorm = res_norm;
128:   KSPLogResidualHistory(ksp,res_norm);
129:   KSPMonitor(ksp,ksp->its,res_norm);

131:   /* check for the convergence - maybe the current guess is good enough */
132:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
133:   if (ksp->reason) {
134:     if (itcount) *itcount = 0;
135:     return(0);
136:   }

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

141:   /* MAIN ITERATION LOOP BEGINNING*/
142:   /* keep iterating until we have converged OR generated the max number
143:      of directions OR reached the max number of iterations for the method */
144:   while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
145:     if (loc_it) {
146:       KSPLogResidualHistory(ksp,res_norm);
147:       KSPMonitor(ksp,ksp->its,res_norm);
148:     }
149:     fgmres->it = (loc_it - 1);

151:     /* see if more space is needed for work vectors */
152:     if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
153:       KSPFGMRESGetNewVectors(ksp,loc_it+1);
154:       /* (loc_it+1) is passed in as number of the first vector that should
155:          be allocated */
156:     }

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


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

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


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

177:     /* new entry in hessenburg is the 2-norm of our new direction */
178:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);

180:     *HH(loc_it+1,loc_it)  = tt;
181:     *HES(loc_it+1,loc_it) = tt;

183:     /* Happy Breakdown Check */
184:     hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
185:     /* RS(loc_it) contains the res_norm from the last iteration  */
186:     hapbnd = PetscMin(fgmres->haptol,hapbnd);
187:     if (tt > hapbnd) {
188:       /* scale new direction by its norm */
189:       VecScale(VEC_VV(loc_it+1),1.0/tt);
190:     } else {
191:       /* This happens when the solution is exactly reached. */
192:       /* So there is no new direction... */
193:       VecSet(VEC_TEMP,0.0);     /* set VEC_TEMP to 0 */
194:       hapend = PETSC_TRUE;
195:     }
196:     /* note that for FGMRES we could get HES(loc_it+1, loc_it)  = 0 and the
197:        current solution would not be exact if HES was singular.  Note that
198:        HH non-singular implies that HES is no singular, and HES is guaranteed
199:        to be nonsingular when PREVECS are linearly independent and A is
200:        nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
201:        of HES). So we should really add a check to verify that HES is nonsingular.*/


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

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

212:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
213:     ksp->its++;
214:     ksp->rnorm = res_norm;
215:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

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

219:     /* Catch error in happy breakdown and signal convergence and break from loop */
220:     if (hapend) {
221:       if (!ksp->reason) {
222:         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_norm);
223:         else {
224:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
225:           break;
226:         }
227:       }
228:     }
229:   }
230:   /* END OF ITERATION LOOP */
231:   KSPLogResidualHistory(ksp,res_norm);

233:   /*
234:      Monitor if we know that we will not return for a restart */
235:   if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) {
236:     KSPMonitor(ksp,ksp->its,res_norm);
237:   }

239:   if (itcount) *itcount = loc_it;

241:   /*
242:     Down here we have to solve for the "best" coefficients of the Krylov
243:     columns, add the solution values together, and possibly unwind the
244:     preconditioning from the solution
245:    */

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

251:   KSPFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
252:   return(0);
253: }

255: /*
256:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

262:    Output Parameter:
263: .     outits - number of iterations used

265: */

267: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
268: {
270:   PetscInt       cycle_its = 0; /* iterations done in a call to KSPFGMRESCycle */
271:   KSP_FGMRES     *fgmres   = (KSP_FGMRES*)ksp->data;
272:   PetscBool      diagonalscale;

275:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
276:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

278:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
279:   ksp->its = 0;
280:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

282:   /* Compute the initial (NOT preconditioned) residual */
283:   if (!ksp->guess_zero) {
284:     KSPFGMRESResidual(ksp);
285:   } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
286:     VecCopy(ksp->vec_rhs,VEC_VV(0));
287:   }
288:   /* now the residual is in VEC_VV(0) - which is what
289:      KSPFGMRESCycle expects... */

291:   KSPFGMRESCycle(&cycle_its,ksp);
292:   while (!ksp->reason) {
293:     KSPFGMRESResidual(ksp);
294:     if (ksp->its >= ksp->max_it) break;
295:     KSPFGMRESCycle(&cycle_its,ksp);
296:   }
297:   /* mark lack of convergence */
298:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
299:   return(0);
300: }

302: extern PetscErrorCode KSPReset_FGMRES(KSP);
303: /*

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

307: */
308: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
309: {

313:   KSPReset_FGMRES(ksp);
314:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",NULL);
315:   KSPDestroy_GMRES(ksp);
316:   return(0);
317: }

319: /*
320:     KSPFGMRESBuildSoln - create the solution from the starting vector and the
321:                       current iterates.

323:     Input parameters:
324:         nrs - work area of size it + 1.
325:         vguess  - index of initial guess
326:         vdest - index of result.  Note that vguess may == vdest (replace
327:                 guess with the solution).
328:         it - HH upper triangular part is a block of size (it+1) x (it+1)

330:      This is an internal routine that knows about the FGMRES internals.
331:  */
332: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
333: {
334:   PetscScalar    tt;
336:   PetscInt       ii,k,j;
337:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)(ksp->data);

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

342:   /* If it is < 0, no fgmres steps have been performed */
343:   if (it < 0) {
344:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
345:     return(0);
346:   }

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

350:   /* solve the upper triangular system - RS is the right side and HH is
351:      the upper triangular matrix  - put soln in nrs */
352:   if (*HH(it,it) != 0.0) {
353:     nrs[it] = *RS(it) / *HH(it,it);
354:   } else {
355:     nrs[it] = 0.0;
356:   }
357:   for (ii=1; ii<=it; ii++) {
358:     k  = it - ii;
359:     tt = *RS(k);
360:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
361:     nrs[k] = tt / *HH(k,k);
362:   }

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

369:   /* put updated solution into vdest.*/
370:   if (vdest != vguess) {
371:     VecCopy(VEC_TEMP,vdest);
372:     VecAXPY(vdest,1.0,vguess);
373:   } else { /* replace guess with solution */
374:     VecAXPY(vdest,1.0,VEC_TEMP);
375:   }
376:   return(0);
377: }

379: /*

381:     KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
382:                             Return new residual.

384:     input parameters:

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

391:     output parameters:
392: .        res - the new residual

394:  */
395: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
396: {
397:   PetscScalar *hh,*cc,*ss,tt;
398:   PetscInt    j;
399:   KSP_FGMRES  *fgmres = (KSP_FGMRES*)(ksp->data);

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

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

412:   for (j=1; j<=it; j++) {
413:     tt  = *hh;
414:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
415:     hh++;
416:     *hh = *cc++ * *hh - (*ss++ * tt);
417:     /* hh, cc, and ss have all been incremented one by end of loop */
418:   }

420:   /*
421:     compute the new plane rotation, and apply it to:
422:      1) the right-hand-side of the Hessenberg system (RS)
423:         note: it affects RS(it) and RS(it+1)
424:      2) the new column of the Hessenberg matrix
425:         note: it affects HH(it,it) which is currently pointed to
426:         by hh and HH(it+1, it) (*(hh+1))
427:     thus obtaining the updated value of the residual...
428:   */

430:   /* compute new plane rotation */

432:   if (!hapend) {
433:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
434:     if (tt == 0.0) {
435:       ksp->reason = KSP_DIVERGED_NULL;
436:       return(0);
437:     }

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

442:     /* apply to 1) and 2) */
443:     *RS(it+1) = -(*ss * *RS(it));
444:     *RS(it)   = PetscConj(*cc) * *RS(it);
445:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);

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

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

457:     *res = 0.0;
458:   }
459:   return(0);
460: }

462: /*

464:    KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
465:                          VEC_VV(it), and more preconditioned work vectors, starting
466:                          from PREVEC(i).

468: */
469: static PetscErrorCode KSPFGMRESGetNewVectors(KSP ksp,PetscInt it)
470: {
471:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
472:   PetscInt       nwork   = fgmres->nwork_alloc; /* number of work vector chunks allocated */
473:   PetscInt       nalloc;                      /* number to allocate */
475:   PetscInt       k;

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

481:   /* Adjust the number to allocate to make sure that we don't exceed the
482:      number of available slots (fgmres->vecs_allocated)*/
483:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated) {
484:     nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
485:   }
486:   if (!nalloc) return(0);

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

490:   /* work vectors */
491:   KSPCreateVecs(ksp,nalloc,&fgmres->user_work[nwork],0,NULL);
492:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
493:   for (k=0; k < nalloc; k++) {
494:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
495:   }
496:   /* specify size of chunk allocated */
497:   fgmres->mwork_alloc[nwork] = nalloc;

499:   /* preconditioned vectors */
500:   KSPCreateVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,NULL);
501:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
502:   for (k=0; k < nalloc; k++) {
503:     fgmres->prevecs[it+k] = fgmres->prevecs_user_work[nwork][k];
504:   }

506:   /* increment the number of work vector chunks */
507:   fgmres->nwork_alloc++;
508:   return(0);
509: }

511: /*

513:    KSPBuildSolution_FGMRES

515:      Input Parameter:
516: .     ksp - the Krylov space object
517: .     ptr-

519:    Output Parameter:
520: .     result - the solution

522:    Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
523:    calls directly.

525: */
526: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
527: {
528:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

532:   if (!ptr) {
533:     if (!fgmres->sol_temp) {
534:       VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
535:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)fgmres->sol_temp);
536:     }
537:     ptr = fgmres->sol_temp;
538:   }
539:   if (!fgmres->nrs) {
540:     /* allocate the work area */
541:     PetscMalloc1(fgmres->max_k,&fgmres->nrs);
542:     PetscLogObjectMemory((PetscObject)ksp,fgmres->max_k*sizeof(PetscScalar));
543:   }

545:   KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
546:   if (result) *result = ptr;
547:   return(0);
548: }

550: PetscErrorCode KSPSetFromOptions_FGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
551: {
553:   PetscBool      flg;

556:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
557:   PetscOptionsHead(PetscOptionsObject,"KSP flexible GMRES Options");
558:   PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
559:   if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
560:   PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
561:   if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
562:   PetscOptionsTail();
563:   return(0);
564: }

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

569: static PetscErrorCode  KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
570: {
573:   ((KSP_FGMRES*)ksp->data)->modifypc      = fcn;
574:   ((KSP_FGMRES*)ksp->data)->modifydestroy = d;
575:   ((KSP_FGMRES*)ksp->data)->modifyctx     = ctx;
576:   return(0);
577: }


580: PetscErrorCode KSPReset_FGMRES(KSP ksp)
581: {
582:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
584:   PetscInt       i;

587:   PetscFree (fgmres->prevecs);
588:   if(fgmres->nwork_alloc>0){
589:     i=0;
590:     /* In the first allocation we allocated VEC_OFFSET fewer vectors in prevecs */
591:     VecDestroyVecs(fgmres->mwork_alloc[i]-VEC_OFFSET,&fgmres->prevecs_user_work[i]);
592:     for (i=1; i<fgmres->nwork_alloc; i++) {
593:       VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
594:     }
595:   }
596:   PetscFree(fgmres->prevecs_user_work);
597:   if (fgmres->modifydestroy) {
598:     (*fgmres->modifydestroy)(fgmres->modifyctx);
599:   }
600:   KSPReset_GMRES(ksp);
601:   return(0);
602: }

604: PetscErrorCode  KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
605: {
606:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;

610:   if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
611:   if (!ksp->setupstage) {
612:     gmres->max_k = max_k;
613:   } else if (gmres->max_k != max_k) {
614:     gmres->max_k    = max_k;
615:     ksp->setupstage = KSP_SETUP_NEW;
616:     /* free the data structures, then create them again */
617:     KSPReset_FGMRES(ksp);
618:   }
619:   return(0);
620: }

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

627:   *max_k = gmres->max_k;
628:   return(0);
629: }

631: /*MC
632:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
633:                 developed by Saad with restart


636:    Options Database Keys:
637: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
638: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
639: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
640:                              vectors are allocated as needed)
641: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
642: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
643: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
644:                                    stability of the classical Gram-Schmidt  orthogonalization.
645: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
646: .   -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
647: -   -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

649:    Level: beginner

651:     Notes:
652:     See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
653:            Only right preconditioning is supported.

655:     Notes:
656:     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)
657:            be bi-CG-stab with a preconditioner of Jacobi.

659:     Developer Notes:
660:     This object is subclassed off of KSPGMRES

662: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
663:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
664:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
665:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
666:            KSPFGMRESModifyPCKSP()

668: M*/

670: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
671: {
672:   KSP_FGMRES     *fgmres;

676:   PetscNewLog(ksp,&fgmres);

678:   ksp->data                              = (void*)fgmres;
679:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
680:   ksp->ops->setup                        = KSPSetUp_FGMRES;
681:   ksp->ops->solve                        = KSPSolve_FGMRES;
682:   ksp->ops->reset                        = KSPReset_FGMRES;
683:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
684:   ksp->ops->view                         = KSPView_GMRES;
685:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
686:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
687:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

689:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
690:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);

692:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
693:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
694:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
695:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_FGMRES);
696:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_FGMRES);
697:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",KSPFGMRESSetModifyPC_FGMRES);
698:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
699:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);


702:   fgmres->haptol         = 1.0e-30;
703:   fgmres->q_preallocate  = 0;
704:   fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
705:   fgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
706:   fgmres->nrs            = 0;
707:   fgmres->sol_temp       = 0;
708:   fgmres->max_k          = FGMRES_DEFAULT_MAXK;
709:   fgmres->Rsvd           = 0;
710:   fgmres->orthogwork     = 0;
711:   fgmres->modifypc       = KSPFGMRESModifyPCNoChange;
712:   fgmres->modifyctx      = NULL;
713:   fgmres->modifydestroy  = NULL;
714:   fgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
715:   return(0);
716: }