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.
 12:  
 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 FGMRESGetNewVectors(KSP,PetscInt);
 21: static PetscErrorCode FGMRESUpdateHessenberg(KSP,PetscInt,PetscBool ,PetscReal *);
 22: static PetscErrorCode BuildFgmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 26: /*

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

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

 33: */
 36: PetscErrorCode    KSPSetUp_FGMRES(KSP ksp)
 37: {
 39:   PetscInt       max_k,k;
 40:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;

 43:   if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"no symmetric preconditioning for KSPFGMRES");
 44:   else if (ksp->pc_side == PC_LEFT) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"no left preconditioning for KSPFGMRES");

 46:   max_k         = fgmres->max_k;

 48:   KSPSetUp_GMRES(ksp);

 50:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs);
 51:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs_user_work);
 52:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)));

 54:   KSPGetVecs(ksp,fgmres->vv_allocated,&fgmres->prevecs_user_work[0],0,PETSC_NULL);
 55:   PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->prevecs_user_work[0]);
 56:   for (k=0; k < fgmres->vv_allocated; k++) {
 57:     fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
 58:   }

 60:   return(0);
 61: }

 63: /* 
 64:     FGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED) 
 65: */
 68: static PetscErrorCode FGMRESResidual(KSP ksp)
 69: {
 70:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);
 71:   Mat            Amat,Pmat;
 72:   MatStructure   pflag;

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

 78:   /* put A*x into VEC_TEMP */
 79:   MatMult(Amat,ksp->vec_sol,VEC_TEMP);
 80:   /* now put residual (-A*x + f) into vec_vv(0) */
 81:   VecWAXPY(VEC_VV(0),-1.0,VEC_TEMP,ksp->vec_rhs);
 82:   return(0);
 83: }

 85: /*

 87:     FGMRESCycle - Run fgmres, possibly with restart.  Return residual 
 88:                   history if requested.

 90:     input parameters:
 91: .         fgmres  - structure containing parameters and work areas

 93:     output parameters:
 94: .        itcount - number of iterations used.  If null, ignored.
 95: .        converged - 0 if not converged

 97:                   
 98:     Notes:
 99:     On entry, the value in vector VEC_VV(0) should be 
100:     the initial residual.


103:  */
106: PetscErrorCode FGMREScycle(PetscInt *itcount,KSP ksp)
107: {

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


121:   /* Number of pseudo iterations since last restart is the number 
122:      of prestart directions */
123:   loc_it = 0;

125:   /* note: (fgmres->it) is always set one less than (loc_it) It is used in 
126:      KSPBUILDSolution_FGMRES, where it is passed to BuildFGmresSoln.  
127:      Note that when BuildFGmresSoln is called from this function, 
128:      (loc_it -1) is passed, so the two are equivalent */
129:   fgmres->it = (loc_it - 1);

131:   /* initial residual is in VEC_VV(0)  - compute its norm*/
132:   VecNorm(VEC_VV(0),NORM_2,&res_norm);

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

138:   ksp->rnorm = res_norm;
139:   KSPLogResidualHistory(ksp,res_norm);

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

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

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

166:     /* CHANGE THE PRECONDITIONER? */
167:     /* ModifyPC is the callback function that can be used to
168:        change the PC or its attributes before its applied */
169:     (*fgmres->modifypc)(ksp,ksp->its,loc_it,res_norm,fgmres->modifyctx);
170: 
171: 
172:     /* apply PRECONDITIONER to direction vector and store with 
173:        preconditioned vectors in prevec */
174:     KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));
175: 
176:     PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
177:     /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
178:     MatMult(Amat,PREVEC(loc_it),VEC_VV(1+loc_it));

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

185:     /* new entry in hessenburg is the 2-norm of our new direction */
186:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
187:     *HH(loc_it+1,loc_it)   = tt;
188:     *HES(loc_it+1,loc_it)  = tt;

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

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

216:     loc_it++;
217:     fgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
218: 
219:     PetscObjectTakeAccess(ksp);
220:     ksp->its++;
221:     ksp->rnorm = res_norm;
222:     PetscObjectGrantAccess(ksp);

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

226:     /* Catch error in happy breakdown and signal convergence and break from loop */
227:     if (hapend) {
228:       if (!ksp->reason) {
229:         SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down,but convergence was not indicated.");
230:       }
231:       break;
232:     }
233:   }
234:   /* END OF ITERATION LOOP */

236:   KSPLogResidualHistory(ksp,res_norm);

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

244:   if (itcount) *itcount    = loc_it;

246:   /*
247:     Down here we have to solve for the "best" coefficients of the Krylov
248:     columns, add the solution values together, and possibly unwind the
249:     preconditioning from the solution
250:    */
251: 
252:   /* Form the solution (or the solution so far) */
253:   /* Note: must pass in (loc_it-1) for iteration count so that BuildFgmresSoln
254:      properly navigates */

256:   BuildFgmresSoln(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 FGMREScycle */
279:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
280:   PetscBool      diagonalscale;

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

286:   PetscObjectTakeAccess(ksp);
287:   ksp->its = 0;
288:   PetscObjectGrantAccess(ksp);

290:   /* Compute the initial (NOT preconditioned) residual */
291:   if (!ksp->guess_zero) {
292:     FGMRESResidual(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:      FGMREScycle expects... */
298: 
299:   FGMREScycle(&cycle_its,ksp);
300:   while (!ksp->reason) {
301:     FGMRESResidual(ksp);
302:     if (ksp->its >= ksp->max_it) break;
303:     FGMREScycle(&cycle_its,ksp);
304:   }
305:   /* mark lack of convergence */
306:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;

308:   return(0);
309: }

313: /*

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

317: */
320: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
321: {

325:   KSPReset_FGMRES(ksp);
326:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C","",PETSC_NULL);
327:   KSPDestroy_GMRES(ksp);
328:   return(0);
329: }

331: /*
332:     BuildFgmresSoln - create the solution from the starting vector and the
333:                       current iterates.

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

342:      This is an internal routine that knows about the FGMRES internals.
343:  */
346: static PetscErrorCode BuildFgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
347: {
348:   PetscScalar    tt;
350:   PetscInt       ii,k,j;
351:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);

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

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

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

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

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

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

393: /*

395:     FGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
396:                             Return new residual.

398:     input parameters:

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

405:     output parameters:
406: .        res - the new residual
407:         
408:  */
411: static PetscErrorCode FGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool  hapend,PetscReal *res)
412: {
413:   PetscScalar   *hh,*cc,*ss,tt;
414:   PetscInt      j;
415:   KSP_FGMRES    *fgmres = (KSP_FGMRES *)(ksp->data);

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

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

428:   for (j=1; j<=it; j++) {
429:     tt  = *hh;
430: #if defined(PETSC_USE_COMPLEX)
431:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
432: #else
433:     *hh = *cc * tt + *ss * *(hh+1);
434: #endif
435:     hh++;
436:     *hh = *cc++ * *hh - (*ss++ * tt);
437:     /* hh, cc, and ss have all been incremented one by end of loop */
438:   }

440:   /*
441:     compute the new plane rotation, and apply it to:
442:      1) the right-hand-side of the Hessenberg system (RS)
443:         note: it affects RS(it) and RS(it+1)
444:      2) the new column of the Hessenberg matrix
445:         note: it affects HH(it,it) which is currently pointed to 
446:         by hh and HH(it+1, it) (*(hh+1))  
447:     thus obtaining the updated value of the residual...
448:   */

450:   /* compute new plane rotation */

452:   if (!hapend) {
453: #if defined(PETSC_USE_COMPLEX)
454:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
455: #else
456:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
457: #endif
458:     if (tt == 0.0) {
459:       ksp->reason = KSP_DIVERGED_NULL;
460:       return(0);
461:     }

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

466:     /* apply to 1) and 2) */
467:     *RS(it+1) = - (*ss * *RS(it));
468: #if defined(PETSC_USE_COMPLEX)
469:     *RS(it)   = PetscConj(*cc) * *RS(it);
470:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);
471: #else
472:     *RS(it)   = *cc * *RS(it);
473:     *hh       = *cc * *hh + *ss * *(hh+1);
474: #endif

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

479:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
480:             another rotation matrix (so RH doesn't change).  The new residual is 
481:             always the new sine term times the residual from last time (RS(it)), 
482:             but now the new sine rotation would be zero...so the residual should
483:             be zero...so we will multiply "zero" by the last residual.  This might
484:             not be exactly what we want to do here -could just return "zero". */
485: 
486:     *res = 0.0;
487:   }
488:   return(0);
489: }

491: /*

493:    FGMRESGetNewVectors - This routine allocates more work vectors, starting from 
494:                          VEC_VV(it), and more preconditioned work vectors, starting 
495:                          from PREVEC(i).

497: */
500: static PetscErrorCode FGMRESGetNewVectors(KSP ksp,PetscInt it)
501: {
502:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
503:   PetscInt       nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
504:   PetscInt       nalloc;                      /* number to allocate */
506:   PetscInt       k;
507: 
509:   nalloc = fgmres->delta_allocate; /* number of vectors to allocate 
510:                                       in a single chunk */

512:   /* Adjust the number to allocate to make sure that we don't exceed the
513:      number of available slots (fgmres->vecs_allocated)*/
514:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated){
515:     nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
516:   }
517:   if (!nalloc) return(0);

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

521:   /* work vectors */
522:   KSPGetVecs(ksp,nalloc,&fgmres->user_work[nwork],0,PETSC_NULL);
523:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
524:   for (k=0; k < nalloc; k++) {
525:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
526:   }
527:   /* specify size of chunk allocated */
528:   fgmres->mwork_alloc[nwork] = nalloc;

530:   /* preconditioned vectors */
531:   KSPGetVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,PETSC_NULL);
532:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
533:   for (k=0; k < nalloc; k++) {
534:     fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
535:   }

537:   /* increment the number of work vector chunks */
538:   fgmres->nwork_alloc++;
539:   return(0);
540: }

542: /* 

544:    KSPBuildSolution_FGMRES

546:      Input Parameter:
547: .     ksp - the Krylov space object
548: .     ptr-

550:    Output Parameter:
551: .     result - the solution

553:    Note: this calls BuildFgmresSoln - the same function that FGMREScycle
554:    calls directly.  

556: */
559: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
560: {
561:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;

565:   if (!ptr) {
566:     if (!fgmres->sol_temp) {
567:       VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
568:       PetscLogObjectParent(ksp,fgmres->sol_temp);
569:     }
570:     ptr = fgmres->sol_temp;
571:   }
572:   if (!fgmres->nrs) {
573:     /* allocate the work area */
574:     PetscMalloc(fgmres->max_k*sizeof(PetscScalar),&fgmres->nrs);
575:     PetscLogObjectMemory(ksp,fgmres->max_k*sizeof(PetscScalar));
576:   }
577: 
578:   BuildFgmresSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
579:   if (result) *result = ptr;
580: 
581:   return(0);
582: }


588: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
589: {
591:   PetscBool      flg;

594:   KSPSetFromOptions_GMRES(ksp);
595:   PetscOptionsHead("KSP flexible GMRES Options");
596:     PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
597:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
598:     PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
599:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
600:   PetscOptionsTail();
601:   return(0);
602: }


608: typedef PetscErrorCode (*FCN2)(void*);
612: PetscErrorCode  KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
613: {
616:   ((KSP_FGMRES *)ksp->data)->modifypc      = fcn;
617:   ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
618:   ((KSP_FGMRES *)ksp->data)->modifyctx     = ctx;
619:   return(0);
620: }



635: PetscErrorCode KSPReset_FGMRES(KSP ksp)
636: {
637:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
639:   PetscInt       i;

642:   PetscFree (fgmres->prevecs);
643:   for (i=0; i<fgmres->nwork_alloc; i++) {
644:     VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
645:   }
646:   PetscFree(fgmres->prevecs_user_work);
647:   if (fgmres->modifydestroy) {
648:     (*fgmres->modifydestroy)(fgmres->modifyctx);
649:   }
650:   KSPReset_GMRES(ksp);
651:   return(0);
652: }

657: PetscErrorCode  KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
658: {
659:   KSP_FGMRES     *gmres = (KSP_FGMRES *)ksp->data;

663:   if (max_k < 1) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
664:   if (!ksp->setupstage) {
665:     gmres->max_k = max_k;
666:   } else if (gmres->max_k != max_k) {
667:      gmres->max_k = max_k;
668:      ksp->setupstage = KSP_SETUP_NEW;
669:      /* free the data structures, then create them again */
670:      KSPReset_FGMRES(ksp);
671:   }
672:   return(0);
673: }

679: PetscErrorCode  KSPGMRESGetRestart_FGMRES(KSP ksp,PetscInt *max_k)
680: {
681:   KSP_FGMRES     *gmres = (KSP_FGMRES *)ksp->data;

684:   *max_k = gmres->max_k;
685:   return(0);
686: }


694: /*MC
695:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.  
696:                 developed by Saad with restart


699:    Options Database Keys:
700: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
701: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
702: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 
703:                              vectors are allocated as needed)
704: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
705: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
706: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 
707:                                    stability of the classical Gram-Schmidt  orthogonalization.
708: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
709: .   -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
710: -   -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

712:    Level: beginner

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

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

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

722: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
723:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
724:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
725:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
726:            KSPFGMRESModifyPCKSP()

728: M*/

733: PetscErrorCode  KSPCreate_FGMRES(KSP ksp)
734: {
735:   KSP_FGMRES     *fgmres;

739:   PetscNewLog(ksp,KSP_FGMRES,&fgmres);
740:   ksp->data                              = (void*)fgmres;
741:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
742:   ksp->ops->setup                        = KSPSetUp_FGMRES;
743:   ksp->ops->solve                        = KSPSolve_FGMRES;
744:   ksp->ops->reset                        = KSPReset_FGMRES;
745:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
746:   ksp->ops->view                         = KSPView_GMRES;
747:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
748:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
749:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

751:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);

753:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
754:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
755:                                      KSPGMRESSetPreAllocateVectors_GMRES);
756:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
757:                                     "KSPGMRESSetOrthogonalization_GMRES",
758:                                      KSPGMRESSetOrthogonalization_GMRES);
759:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
760:                                     "KSPGMRESGetOrthogonalization_GMRES",
761:                                      KSPGMRESGetOrthogonalization_GMRES);
762:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
763:                                     "KSPGMRESSetRestart_FGMRES",
764:                                      KSPGMRESSetRestart_FGMRES);
765:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
766:                                     "KSPGMRESGetRestart_FGMRES",
767:                                      KSPGMRESGetRestart_FGMRES);
768:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
769:                                     "KSPFGMRESSetModifyPC_FGMRES",
770:                                      KSPFGMRESSetModifyPC_FGMRES);
771:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
772:                                     "KSPGMRESSetCGSRefinementType_GMRES",
773:                                      KSPGMRESSetCGSRefinementType_GMRES);
774:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
775:                                     "KSPGMRESGetCGSRefinementType_GMRES",
776:                                      KSPGMRESGetCGSRefinementType_GMRES);


779:   fgmres->haptol              = 1.0e-30;
780:   fgmres->q_preallocate       = 0;
781:   fgmres->delta_allocate      = FGMRES_DELTA_DIRECTIONS;
782:   fgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
783:   fgmres->nrs                 = 0;
784:   fgmres->sol_temp            = 0;
785:   fgmres->max_k               = FGMRES_DEFAULT_MAXK;
786:   fgmres->Rsvd                = 0;
787:   fgmres->orthogwork          = 0;
788:   fgmres->modifypc            = KSPFGMRESModifyPCNoChange;
789:   fgmres->modifyctx           = PETSC_NULL;
790:   fgmres->modifydestroy       = PETSC_NULL;
791:   fgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
792:   return(0);
793: }