Actual source code: fgmres.c

petsc-3.7.7 2017-09-25
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>       /*I  "petscksp.h"  I*/
 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: */
 34: PetscErrorCode    KSPSetUp_FGMRES(KSP ksp)
 35: {
 37:   PetscInt       max_k,k;
 38:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

 41:   max_k = fgmres->max_k;

 43:   KSPSetUp_GMRES(ksp);

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

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

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

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

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

 81: /*

 83:     KSPFGMRESCycle - Run fgmres, possibly with restart.  Return residual
 84:                   history if requested.

 86:     input parameters:
 87: .        fgmres  - structure containing parameters and work areas

 89:     output parameters:
 90: .        itcount - number of iterations used.  If null, ignored.
 91: .        converged - 0 if not converged


 94:     Notes:
 95:     On entry, the value in vector VEC_VV(0) should be
 96:     the initial residual.


 99:  */
102: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount,KSP ksp)
103: {

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

115:   /* Number of pseudo iterations since last restart is the number
116:      of prestart directions */
117:   loc_it = 0;

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

125:   /* initial residual is in VEC_VV(0)  - compute its norm*/
126:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
127:   KSPCheckNorm(ksp,res_norm);

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

133:   ksp->rnorm = res_norm;
134:   KSPLogResidualHistory(ksp,res_norm);
135:   KSPMonitor(ksp,ksp->its,res_norm);

137:   /* check for the convergence - maybe the current guess is good enough */
138:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
139:   if (ksp->reason) {
140:     if (itcount) *itcount = 0;
141:     return(0);
142:   }

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

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

157:     /* see if more space is needed for work vectors */
158:     if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
159:       KSPFGMRESGetNewVectors(ksp,loc_it+1);
160:       /* (loc_it+1) is passed in as number of the first vector that should
161:          be allocated */
162:     }

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


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

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


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

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

186:     *HH(loc_it+1,loc_it)  = tt;
187:     *HES(loc_it+1,loc_it) = tt;

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


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

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

218:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
219:     ksp->its++;
220:     ksp->rnorm = res_norm;
221:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

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

225:     /* Catch error in happy breakdown and signal convergence and break from loop */
226:     if (hapend) {
227:       if (!ksp->reason) {
228:         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);
229:         else {
230:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
231:           break;
232:         }
233:       }
234:     }
235:   }
236:   /* END OF ITERATION LOOP */
237:   KSPLogResidualHistory(ksp,res_norm);

239:   /*
240:      Monitor if we know that we will not return for a restart */
241:   if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) {
242:     KSPMonitor(ksp,ksp->its,res_norm);
243:   }

245:   if (itcount) *itcount = loc_it;

247:   /*
248:     Down here we have to solve for the "best" coefficients of the Krylov
249:     columns, add the solution values together, and possibly unwind the
250:     preconditioning from the solution
251:    */

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

257:   KSPFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
258:   return(0);
259: }

261: /*
262:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

268:    Output Parameter:
269: .     outits - number of iterations used

271: */

275: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
276: {
278:   PetscInt       cycle_its = 0; /* iterations done in a call to KSPFGMRESCycle */
279:   KSP_FGMRES     *fgmres   = (KSP_FGMRES*)ksp->data;
280:   PetscBool      diagonalscale;

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

286:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
287:   ksp->its = 0;
288:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

290:   /* Compute the initial (NOT preconditioned) residual */
291:   if (!ksp->guess_zero) {
292:     KSPFGMRESResidual(ksp);
293:   } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
294:     VecCopy(ksp->vec_rhs,VEC_VV(0));
295:   }
296:   /* now the residual is in VEC_VV(0) - which is what
297:      KSPFGMRESCycle expects... */

299:   KSPFGMRESCycle(&cycle_its,ksp);
300:   while (!ksp->reason) {
301:     KSPFGMRESResidual(ksp);
302:     if (ksp->its >= ksp->max_it) break;
303:     KSPFGMRESCycle(&cycle_its,ksp);
304:   }
305:   /* mark lack of convergence */
306:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
307:   return(0);
308: }

310: extern PetscErrorCode KSPReset_FGMRES(KSP);
311: /*

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

315: */
318: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
319: {

323:   KSPReset_FGMRES(ksp);
324:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",NULL);
325:   KSPDestroy_GMRES(ksp);
326:   return(0);
327: }

329: /*
330:     KSPFGMRESBuildSoln - create the solution from the starting vector and the
331:                       current iterates.

333:     Input parameters:
334:         nrs - work area of size it + 1.
335:         vguess  - index of initial guess
336:         vdest - index of result.  Note that vguess may == vdest (replace
337:                 guess with the solution).
338:         it - HH upper triangular part is a block of size (it+1) x (it+1)

340:      This is an internal routine that knows about the FGMRES internals.
341:  */
344: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
345: {
346:   PetscScalar    tt;
348:   PetscInt       ii,k,j;
349:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)(ksp->data);

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

354:   /* If it is < 0, no fgmres steps have been performed */
355:   if (it < 0) {
356:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
357:     return(0);
358:   }

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

362:   /* solve the upper triangular system - RS is the right side and HH is
363:      the upper triangular matrix  - put soln in nrs */
364:   if (*HH(it,it) != 0.0) {
365:     nrs[it] = *RS(it) / *HH(it,it);
366:   } else {
367:     nrs[it] = 0.0;
368:   }
369:   for (ii=1; ii<=it; ii++) {
370:     k  = it - ii;
371:     tt = *RS(k);
372:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
373:     nrs[k] = tt / *HH(k,k);
374:   }

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

381:   /* put updated solution into vdest.*/
382:   if (vdest != vguess) {
383:     VecCopy(VEC_TEMP,vdest);
384:     VecAXPY(vdest,1.0,vguess);
385:   } else { /* replace guess with solution */
386:     VecAXPY(vdest,1.0,VEC_TEMP);
387:   }
388:   return(0);
389: }

391: /*

393:     KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
394:                             Return new residual.

396:     input parameters:

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

403:     output parameters:
404: .        res - the new residual

406:  */
409: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
410: {
411:   PetscScalar *hh,*cc,*ss,tt;
412:   PetscInt    j;
413:   KSP_FGMRES  *fgmres = (KSP_FGMRES*)(ksp->data);

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

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

426:   for (j=1; j<=it; j++) {
427:     tt  = *hh;
428:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
429:     hh++;
430:     *hh = *cc++ * *hh - (*ss++ * tt);
431:     /* hh, cc, and ss have all been incremented one by end of loop */
432:   }

434:   /*
435:     compute the new plane rotation, and apply it to:
436:      1) the right-hand-side of the Hessenberg system (RS)
437:         note: it affects RS(it) and RS(it+1)
438:      2) the new column of the Hessenberg matrix
439:         note: it affects HH(it,it) which is currently pointed to
440:         by hh and HH(it+1, it) (*(hh+1))
441:     thus obtaining the updated value of the residual...
442:   */

444:   /* compute new plane rotation */

446:   if (!hapend) {
447:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
448:     if (tt == 0.0) {
449:       ksp->reason = KSP_DIVERGED_NULL;
450:       return(0);
451:     }

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

456:     /* apply to 1) and 2) */
457:     *RS(it+1) = -(*ss * *RS(it));
458:     *RS(it)   = PetscConj(*cc) * *RS(it);
459:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);

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

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

471:     *res = 0.0;
472:   }
473:   return(0);
474: }

476: /*

478:    KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
479:                          VEC_VV(it), and more preconditioned work vectors, starting
480:                          from PREVEC(i).

482: */
485: static PetscErrorCode KSPFGMRESGetNewVectors(KSP ksp,PetscInt it)
486: {
487:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
488:   PetscInt       nwork   = fgmres->nwork_alloc; /* number of work vector chunks allocated */
489:   PetscInt       nalloc;                      /* number to allocate */
491:   PetscInt       k;

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

497:   /* Adjust the number to allocate to make sure that we don't exceed the
498:      number of available slots (fgmres->vecs_allocated)*/
499:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated) {
500:     nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
501:   }
502:   if (!nalloc) return(0);

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

506:   /* work vectors */
507:   KSPCreateVecs(ksp,nalloc,&fgmres->user_work[nwork],0,NULL);
508:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
509:   for (k=0; k < nalloc; k++) {
510:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
511:   }
512:   /* specify size of chunk allocated */
513:   fgmres->mwork_alloc[nwork] = nalloc;

515:   /* preconditioned vectors */
516:   KSPCreateVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,NULL);
517:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
518:   for (k=0; k < nalloc; k++) {
519:     fgmres->prevecs[it+k] = fgmres->prevecs_user_work[nwork][k];
520:   }

522:   /* increment the number of work vector chunks */
523:   fgmres->nwork_alloc++;
524:   return(0);
525: }

527: /*

529:    KSPBuildSolution_FGMRES

531:      Input Parameter:
532: .     ksp - the Krylov space object
533: .     ptr-

535:    Output Parameter:
536: .     result - the solution

538:    Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
539:    calls directly.

541: */
544: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
545: {
546:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

550:   if (!ptr) {
551:     if (!fgmres->sol_temp) {
552:       VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
553:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)fgmres->sol_temp);
554:     }
555:     ptr = fgmres->sol_temp;
556:   }
557:   if (!fgmres->nrs) {
558:     /* allocate the work area */
559:     PetscMalloc1(fgmres->max_k,&fgmres->nrs);
560:     PetscLogObjectMemory((PetscObject)ksp,fgmres->max_k*sizeof(PetscScalar));
561:   }

563:   KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
564:   if (result) *result = ptr;
565:   return(0);
566: }

570: PetscErrorCode KSPSetFromOptions_FGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
571: {
573:   PetscBool      flg;

576:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
577:   PetscOptionsHead(PetscOptionsObject,"KSP flexible GMRES Options");
578:   PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
579:   if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
580:   PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
581:   if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
582:   PetscOptionsTail();
583:   return(0);
584: }

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

591: static PetscErrorCode  KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
592: {
595:   ((KSP_FGMRES*)ksp->data)->modifypc      = fcn;
596:   ((KSP_FGMRES*)ksp->data)->modifydestroy = d;
597:   ((KSP_FGMRES*)ksp->data)->modifyctx     = ctx;
598:   return(0);
599: }


604: PetscErrorCode KSPReset_FGMRES(KSP ksp)
605: {
606:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
608:   PetscInt       i;

611:   PetscFree (fgmres->prevecs);
612:   if(fgmres->nwork_alloc>0){
613:     i=0;
614:     /* In the first allocation we allocated VEC_OFFSET fewer vectors in prevecs */
615:     VecDestroyVecs(fgmres->mwork_alloc[i]-VEC_OFFSET,&fgmres->prevecs_user_work[i]);
616:     for (i=1; i<fgmres->nwork_alloc; i++) {
617:       VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
618:     }
619:   }
620:   PetscFree(fgmres->prevecs_user_work);
621:   if (fgmres->modifydestroy) {
622:     (*fgmres->modifydestroy)(fgmres->modifyctx);
623:   }
624:   KSPReset_GMRES(ksp);
625:   return(0);
626: }

630: PetscErrorCode  KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
631: {
632:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;

636:   if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
637:   if (!ksp->setupstage) {
638:     gmres->max_k = max_k;
639:   } else if (gmres->max_k != max_k) {
640:     gmres->max_k    = max_k;
641:     ksp->setupstage = KSP_SETUP_NEW;
642:     /* free the data structures, then create them again */
643:     KSPReset_FGMRES(ksp);
644:   }
645:   return(0);
646: }

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

655:   *max_k = gmres->max_k;
656:   return(0);
657: }

659: /*MC
660:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
661:                 developed by Saad with restart


664:    Options Database Keys:
665: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
666: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
667: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
668:                              vectors are allocated as needed)
669: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
670: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
671: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
672:                                    stability of the classical Gram-Schmidt  orthogonalization.
673: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
674: .   -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
675: -   -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

677:    Level: beginner

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

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

685:     Developer Notes: This object is subclassed off of KSPGMRES

687: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
688:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
689:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
690:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
691:            KSPFGMRESModifyPCKSP()

693: M*/

697: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
698: {
699:   KSP_FGMRES     *fgmres;

703:   PetscNewLog(ksp,&fgmres);

705:   ksp->data                              = (void*)fgmres;
706:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
707:   ksp->ops->setup                        = KSPSetUp_FGMRES;
708:   ksp->ops->solve                        = KSPSolve_FGMRES;
709:   ksp->ops->reset                        = KSPReset_FGMRES;
710:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
711:   ksp->ops->view                         = KSPView_GMRES;
712:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
713:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
714:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

716:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
717:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,0);

719:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
720:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
721:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
722:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_FGMRES);
723:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_FGMRES);
724:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",KSPFGMRESSetModifyPC_FGMRES);
725:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
726:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);


729:   fgmres->haptol         = 1.0e-30;
730:   fgmres->q_preallocate  = 0;
731:   fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
732:   fgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
733:   fgmres->nrs            = 0;
734:   fgmres->sol_temp       = 0;
735:   fgmres->max_k          = FGMRES_DEFAULT_MAXK;
736:   fgmres->Rsvd           = 0;
737:   fgmres->orthogwork     = 0;
738:   fgmres->modifypc       = KSPFGMRESModifyPCNoChange;
739:   fgmres->modifyctx      = NULL;
740:   fgmres->modifydestroy  = NULL;
741:   fgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
742:   return(0);
743: }