Actual source code: fgmres.c

  1: #define PETSCKSP_DLL

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

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

 12:     Restarts:  Restarts are basically solves with x0 not equal to zero.
 13:  
 14:        Contributed by Allison Baker

 16: */

 18:  #include ../src/ksp/ksp/impls/gmres/fgmres/fgmresimpl.h
 19: #define FGMRES_DELTA_DIRECTIONS 10
 20: #define FGMRES_DEFAULT_MAXK     30
 21: static PetscErrorCode FGMRESGetNewVectors(KSP,PetscInt);
 22: static PetscErrorCode FGMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal *);
 23: static PetscErrorCode BuildFgmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 25: EXTERN PetscErrorCode KSPView_GMRES(KSP,PetscViewer);
 26: EXTERN PetscErrorCode KSPSetUp_GMRES(KSP);
 27: /*

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

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

 34: */
 37: PetscErrorCode    KSPSetUp_FGMRES(KSP ksp)
 38: {
 40:   PetscInt       max_k,k;
 41:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;

 44:   if (ksp->pc_side == PC_SYMMETRIC) {
 45:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPFGMRES");
 46:   } else if (ksp->pc_side == PC_LEFT) {
 47:     SETERRQ(PETSC_ERR_SUP,"no left preconditioning for KSPFGMRES");
 48:   }
 49:   max_k         = fgmres->max_k;

 51:   KSPSetUp_GMRES(ksp);

 53:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs);
 54:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs_user_work);
 55:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)));

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

 63:   return(0);
 64: }

 66: /* 
 67:     FGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED) 
 68: */
 71: static PetscErrorCode FGMRESResidual(KSP ksp)
 72: {
 73:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);
 74:   Mat            Amat,Pmat;
 75:   MatStructure   pflag;

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

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

 88: /*

 90:     FGMRESCycle - Run fgmres, possibly with restart.  Return residual 
 91:                   history if requested.

 93:     input parameters:
 94: .         fgmres  - structure containing parameters and work areas

 96:     output parameters:
 97: .        itcount - number of iterations used.  If null, ignored.
 98: .        converged - 0 if not converged

100:                   
101:     Notes:
102:     On entry, the value in vector VEC_VV(0) should be 
103:     the initial residual.


106:  */
109: PetscErrorCode FGMREScycle(PetscInt *itcount,KSP ksp)
110: {

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


124:   /* Number of pseudo iterations since last restart is the number 
125:      of prestart directions */
126:   loc_it = 0;

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

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

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

141:   ksp->rnorm = res_norm;
142:   KSPLogResidualHistory(ksp,res_norm);

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

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

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

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

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

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

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

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

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

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

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

239:   KSPLogResidualHistory(ksp,res_norm);

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

247:   if (itcount) *itcount    = loc_it;

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

259:   BuildFgmresSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);

261:   return(0);
262: }

264: /*  
265:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

271:    Output Parameter:
272: .     outits - number of iterations used

274: */

278: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
279: {
281:   PetscInt       cycle_its = 0; /* iterations done in a call to FGMREScycle */
282:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
283:   PetscTruth     diagonalscale;

286:   PCDiagonalScale(ksp->pc,&diagonalscale);
287:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
288:   if (ksp->normtype != KSP_NORM_UNPRECONDITIONED) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Can only use FGMRES with unpreconditioned residual (it is coded with right preconditioning)");

290:   PetscObjectTakeAccess(ksp);
291:   ksp->its = 0;
292:   PetscObjectGrantAccess(ksp);

294:   /* Compute the initial (NOT preconditioned) residual */
295:   if (!ksp->guess_zero) {
296:     FGMRESResidual(ksp);
297:   } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
298:     VecCopy(ksp->vec_rhs,VEC_VV(0));
299:   }
300:   /* now the residual is in VEC_VV(0) - which is what 
301:      FGMREScycle expects... */
302: 
303:   FGMREScycle(&cycle_its,ksp);
304:   while (!ksp->reason) {
305:     FGMRESResidual(ksp);
306:     if (ksp->its >= ksp->max_it) break;
307:     FGMREScycle(&cycle_its,ksp);
308:   }
309:   /* mark lack of convergence */
310:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;

312:   return(0);
313: }

316: /*

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

320: */
323: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
324: {
325:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
327:   PetscInt       i;

330:   PetscFree (fgmres->prevecs);
331:   for (i=0; i<fgmres->nwork_alloc; i++) {
332:     VecDestroyVecs(fgmres->prevecs_user_work[i],fgmres->mwork_alloc[i]);
333:   }
334:   PetscFree(fgmres->prevecs_user_work);
335:   if (fgmres->modifydestroy) {
336:     (*fgmres->modifydestroy)(fgmres->modifyctx);
337:   }

339:   /* clear composed functions */
340:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C","",PETSC_NULL);
341:   KSPDestroy_GMRES(ksp);
342:   return(0);
343: }

345: /*
346:     BuildFgmresSoln - create the solution from the starting vector and the
347:                       current iterates.

349:     Input parameters:
350:         nrs - work area of size it + 1.
351:         vguess  - index of initial guess
352:         vdest - index of result.  Note that vguess may == vdest (replace
353:                 guess with the solution).
354:         it - HH upper triangular part is a block of size (it+1) x (it+1)  

356:      This is an internal routine that knows about the FGMRES internals.
357:  */
360: static PetscErrorCode BuildFgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
361: {
362:   PetscScalar    tt;
364:   PetscInt       ii,k,j;
365:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);

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

370:   /* If it is < 0, no fgmres steps have been performed */
371:   if (it < 0) {
372:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
373:     return(0);
374:   }

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

378:   /* solve the upper triangular system - RS is the right side and HH is 
379:      the upper triangular matrix  - put soln in nrs */
380:   if (*HH(it,it) != 0.0) {
381:     nrs[it] = *RS(it) / *HH(it,it);
382:   } else {
383:     nrs[it] = 0.0;
384:   }
385:   for (ii=1; ii<=it; ii++) {
386:     k   = it - ii;
387:     tt  = *RS(k);
388:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
389:     nrs[k]   = tt / *HH(k,k);
390:   }

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

397:   /* put updated solution into vdest.*/
398:   if (vdest != vguess) {
399:     VecCopy(VEC_TEMP,vdest);
400:     VecAXPY(vdest,1.0,vguess);
401:   } else  {/* replace guess with solution */
402:     VecAXPY(vdest,1.0,VEC_TEMP);
403:   }
404:   return(0);
405: }

407: /*

409:     FGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
410:                             Return new residual.

412:     input parameters:

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

419:     output parameters:
420: .        res - the new residual
421:         
422:  */
425: static PetscErrorCode FGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
426: {
427:   PetscScalar   *hh,*cc,*ss,tt;
428:   PetscInt      j;
429:   KSP_FGMRES    *fgmres = (KSP_FGMRES *)(ksp->data);

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

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

442:   for (j=1; j<=it; j++) {
443:     tt  = *hh;
444: #if defined(PETSC_USE_COMPLEX)
445:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
446: #else
447:     *hh = *cc * tt + *ss * *(hh+1);
448: #endif
449:     hh++;
450:     *hh = *cc++ * *hh - (*ss++ * tt);
451:     /* hh, cc, and ss have all been incremented one by end of loop */
452:   }

454:   /*
455:     compute the new plane rotation, and apply it to:
456:      1) the right-hand-side of the Hessenberg system (RS)
457:         note: it affects RS(it) and RS(it+1)
458:      2) the new column of the Hessenberg matrix
459:         note: it affects HH(it,it) which is currently pointed to 
460:         by hh and HH(it+1, it) (*(hh+1))  
461:     thus obtaining the updated value of the residual...
462:   */

464:   /* compute new plane rotation */

466:   if (!hapend) {
467: #if defined(PETSC_USE_COMPLEX)
468:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
469: #else
470:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
471: #endif
472:     if (tt == 0.0) {
473:       ksp->reason = KSP_DIVERGED_NULL;
474:       return(0);
475:     }

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

480:     /* apply to 1) and 2) */
481:     *RS(it+1) = - (*ss * *RS(it));
482: #if defined(PETSC_USE_COMPLEX)
483:     *RS(it)   = PetscConj(*cc) * *RS(it);
484:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);
485: #else
486:     *RS(it)   = *cc * *RS(it);
487:     *hh       = *cc * *hh + *ss * *(hh+1);
488: #endif

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

493:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
494:             another rotation matrix (so RH doesn't change).  The new residual is 
495:             always the new sine term times the residual from last time (RS(it)), 
496:             but now the new sine rotation would be zero...so the residual should
497:             be zero...so we will multiply "zero" by the last residual.  This might
498:             not be exactly what we want to do here -could just return "zero". */
499: 
500:     *res = 0.0;
501:   }
502:   return(0);
503: }

505: /*

507:    FGMRESGetNewVectors - This routine allocates more work vectors, starting from 
508:                          VEC_VV(it), and more preconditioned work vectors, starting 
509:                          from PREVEC(i).

511: */
514: static PetscErrorCode FGMRESGetNewVectors(KSP ksp,PetscInt it)
515: {
516:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
517:   PetscInt       nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
518:   PetscInt       nalloc;                      /* number to allocate */
520:   PetscInt       k;
521: 
523:   nalloc = fgmres->delta_allocate; /* number of vectors to allocate 
524:                                       in a single chunk */

526:   /* Adjust the number to allocate to make sure that we don't exceed the
527:      number of available slots (fgmres->vecs_allocated)*/
528:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated){
529:     nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
530:   }
531:   if (!nalloc) return(0);

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

535:   /* work vectors */
536:   KSPGetVecs(ksp,nalloc,&fgmres->user_work[nwork],0,PETSC_NULL);
537:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
538:   for (k=0; k < nalloc; k++) {
539:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
540:   }
541:   /* specify size of chunk allocated */
542:   fgmres->mwork_alloc[nwork] = nalloc;

544:   /* preconditioned vectors */
545:   KSPGetVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,PETSC_NULL);
546:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
547:   for (k=0; k < nalloc; k++) {
548:     fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
549:   }

551:   /* increment the number of work vector chunks */
552:   fgmres->nwork_alloc++;
553:   return(0);
554: }

556: /* 

558:    KSPBuildSolution_FGMRES

560:      Input Parameter:
561: .     ksp - the Krylov space object
562: .     ptr-

564:    Output Parameter:
565: .     result - the solution

567:    Note: this calls BuildFgmresSoln - the same function that FGMREScycle
568:    calls directly.  

570: */
573: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
574: {
575:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;

579:   if (!ptr) {
580:     if (!fgmres->sol_temp) {
581:       VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
582:       PetscLogObjectParent(ksp,fgmres->sol_temp);
583:     }
584:     ptr = fgmres->sol_temp;
585:   }
586:   if (!fgmres->nrs) {
587:     /* allocate the work area */
588:     PetscMalloc(fgmres->max_k*sizeof(PetscScalar),&fgmres->nrs);
589:     PetscLogObjectMemory(ksp,fgmres->max_k*sizeof(PetscScalar));
590:   }
591: 
592:   BuildFgmresSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
593:   if (result) *result = ptr;
594: 
595:   return(0);
596: }


602: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
603: {
605:   PetscTruth     flg;

608:   KSPSetFromOptions_GMRES(ksp);
609:   PetscOptionsHead("KSP flexible GMRES Options");
610:     PetscOptionsTruthGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
611:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
612:     PetscOptionsTruthGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
613:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
614:   PetscOptionsTail();
615:   return(0);
616: }

618: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
619: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);

622: typedef PetscErrorCode (*FCN2)(void*);
626: PetscErrorCode  KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
627: {
630:   ((KSP_FGMRES *)ksp->data)->modifypc      = fcn;
631:   ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
632:   ((KSP_FGMRES *)ksp->data)->modifyctx     = ctx;
633:   return(0);
634: }

638: EXTERN PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP);
639: EXTERN PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP,PetscInt);
640: EXTERN PetscErrorCode  KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));

643: EXTERN PetscErrorCode KSPDestroy_GMRES_Internal(KSP);

647: PetscErrorCode KSPDestroy_FGMRES_Internal(KSP ksp)
648: {
649:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;

653:   KSPDestroy_GMRES_Internal(ksp);
654:   PetscFree (gmres->prevecs);
655:   PetscFree(gmres->prevecs_user_work);
656:   if (gmres->modifydestroy) {
657:     (*gmres->modifydestroy)(gmres->modifyctx);
658:   }
659:   return(0);
660: }

665: PetscErrorCode  KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
666: {
667:   KSP_FGMRES     *gmres = (KSP_FGMRES *)ksp->data;

671:   if (max_k < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
672:   if (!ksp->setupcalled) {
673:     gmres->max_k = max_k;
674:   } else if (gmres->max_k != max_k) {
675:      gmres->max_k = max_k;
676:      ksp->setupcalled = 0;
677:      /* free the data structures, then create them again */
678:      KSPDestroy_FGMRES_Internal(ksp);
679:   }
680:   return(0);
681: }

685: EXTERN PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);

688: /*MC
689:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.  
690:                 developed by Saad with restart


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

706:    Level: beginner

708:     Notes: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
709:            This object is subclassed off of KSPGMRES

711: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
712:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
713:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
714:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
715:            KSPFGMRESModifyPCKSP()

717: M*/

722: PetscErrorCode  KSPCreate_FGMRES(KSP ksp)
723: {
724:   KSP_FGMRES     *fgmres;

728:   PetscNewLog(ksp,KSP_FGMRES,&fgmres);
729:   ksp->data                              = (void*)fgmres;
730:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
731:   ksp->ops->setup                        = KSPSetUp_FGMRES;
732:   ksp->ops->solve                        = KSPSolve_FGMRES;
733:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
734:   ksp->ops->view                         = KSPView_GMRES;
735:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
736:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
737:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

739:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
740:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
741:                                      KSPGMRESSetPreAllocateVectors_GMRES);
742:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
743:                                     "KSPGMRESSetOrthogonalization_GMRES",
744:                                      KSPGMRESSetOrthogonalization_GMRES);
745:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
746:                                     "KSPGMRESSetRestart_FGMRES",
747:                                      KSPGMRESSetRestart_FGMRES);
748:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
749:                                     "KSPFGMRESSetModifyPC_FGMRES",
750:                                      KSPFGMRESSetModifyPC_FGMRES);
751:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
752:                                     "KSPGMRESSetCGSRefinementType_GMRES",
753:                                      KSPGMRESSetCGSRefinementType_GMRES);


756:   fgmres->haptol              = 1.0e-30;
757:   fgmres->q_preallocate       = 0;
758:   fgmres->delta_allocate      = FGMRES_DELTA_DIRECTIONS;
759:   fgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
760:   fgmres->nrs                 = 0;
761:   fgmres->sol_temp            = 0;
762:   fgmres->max_k               = FGMRES_DEFAULT_MAXK;
763:   fgmres->Rsvd                = 0;
764:   fgmres->orthogwork          = 0;
765:   fgmres->modifypc            = KSPFGMRESModifyPCNoChange;
766:   fgmres->modifyctx           = PETSC_NULL;
767:   fgmres->modifydestroy       = PETSC_NULL;
768:   fgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
769:   /*
770:         This is not great since it changes this without explicit request from the user
771:      but there is no left preconditioning in the FGMRES
772:   */
773:   PetscInfo(ksp,"WARNING! Setting PC_SIDE for FGMRES to right!\n");
774:   ksp->pc_side  = PC_RIGHT;
775:   ksp->normtype = KSP_NORM_UNPRECONDITIONED;
776:   return(0);
777: }