Actual source code: fgmres.c

petsc-3.4.4 2014-03-13
  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:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs);
 46:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs_user_work);
 47:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)));

 49:   KSPGetVecs(ksp,fgmres->vv_allocated,&fgmres->prevecs_user_work[0],0,NULL);
 50:   PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->prevecs_user_work[0]);
 51:   for (k=0; k < fgmres->vv_allocated; k++) {
 52:     fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
 53:   }
 54:   return(0);
 55: }

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

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

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

 79: /*

 81:     KSPFGMRESCycle - Run fgmres, possibly with restart.  Return residual
 82:                   history if requested.

 84:     input parameters:
 85: .        fgmres  - structure containing parameters and work areas

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


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


 97:  */
100: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount,KSP ksp)
101: {

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

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

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

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

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

131:   ksp->rnorm = res_norm;
132:   KSPLogResidualHistory(ksp,res_norm);
133:   KSPMonitor(ksp,ksp->its,res_norm);

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

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

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

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

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


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

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


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

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

184:     *HH(loc_it+1,loc_it)  = tt;
185:     *HES(loc_it+1,loc_it) = tt;

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


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

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

216:     PetscObjectAMSTakeAccess((PetscObject)ksp);
217:     ksp->its++;
218:     ksp->rnorm = res_norm;
219:     PetscObjectAMSGrantAccess((PetscObject)ksp);

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

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

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

243:   if (itcount) *itcount = loc_it;

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

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

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

259: /*
260:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

266:    Output Parameter:
267: .     outits - number of iterations used

269: */

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

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

284:   PetscObjectAMSTakeAccess((PetscObject)ksp);
285:   ksp->its = 0;
286:   PetscObjectAMSGrantAccess((PetscObject)ksp);

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

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

308: extern PetscErrorCode KSPReset_FGMRES(KSP);
309: /*

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

313: */
316: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
317: {

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

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

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

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

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

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

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

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

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

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

389: /*

391:     KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
392:                             Return new residual.

394:     input parameters:

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

401:     output parameters:
402: .        res - the new residual

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

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

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

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

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

442:   /* compute new plane rotation */

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

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

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

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

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

469:     *res = 0.0;
470:   }
471:   return(0);
472: }

474: /*

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

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

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

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

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

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

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

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

525: /*

527:    KSPBuildSolution_FGMRES

529:      Input Parameter:
530: .     ksp - the Krylov space object
531: .     ptr-

533:    Output Parameter:
534: .     result - the solution

536:    Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
537:    calls directly.

539: */
542: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
543: {
544:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

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

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

568: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
569: {
571:   PetscBool      flg;

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

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

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


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

609:   PetscFree (fgmres->prevecs);
610:   for (i=0; i<fgmres->nwork_alloc; i++) {
611:     VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
612:   }
613:   PetscFree(fgmres->prevecs_user_work);
614:   if (fgmres->modifydestroy) {
615:     (*fgmres->modifydestroy)(fgmres->modifyctx);
616:   }
617:   KSPReset_GMRES(ksp);
618:   return(0);
619: }

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

629:   if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
630:   if (!ksp->setupstage) {
631:     gmres->max_k = max_k;
632:   } else if (gmres->max_k != max_k) {
633:     gmres->max_k    = max_k;
634:     ksp->setupstage = KSP_SETUP_NEW;
635:     /* free the data structures, then create them again */
636:     KSPReset_FGMRES(ksp);
637:   }
638:   return(0);
639: }

643: PetscErrorCode  KSPGMRESGetRestart_FGMRES(KSP ksp,PetscInt *max_k)
644: {
645:   KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;

648:   *max_k = gmres->max_k;
649:   return(0);
650: }

652: /*MC
653:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
654:                 developed by Saad with restart


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

670:    Level: beginner

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

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

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

680: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
681:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
682:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
683:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
684:            KSPFGMRESModifyPCKSP()

686: M*/

690: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
691: {
692:   KSP_FGMRES     *fgmres;

696:   PetscNewLog(ksp,KSP_FGMRES,&fgmres);

698:   ksp->data                              = (void*)fgmres;
699:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;
700:   ksp->ops->setup                        = KSPSetUp_FGMRES;
701:   ksp->ops->solve                        = KSPSolve_FGMRES;
702:   ksp->ops->reset                        = KSPReset_FGMRES;
703:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
704:   ksp->ops->view                         = KSPView_GMRES;
705:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
706:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
707:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

709:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
710:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,0);

712:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
713:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
714:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
715:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_FGMRES);
716:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_FGMRES);
717:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",KSPFGMRESSetModifyPC_FGMRES);
718:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
719:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);


722:   fgmres->haptol         = 1.0e-30;
723:   fgmres->q_preallocate  = 0;
724:   fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
725:   fgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
726:   fgmres->nrs            = 0;
727:   fgmres->sol_temp       = 0;
728:   fgmres->max_k          = FGMRES_DEFAULT_MAXK;
729:   fgmres->Rsvd           = 0;
730:   fgmres->orthogwork     = 0;
731:   fgmres->modifypc       = KSPFGMRESModifyPCNoChange;
732:   fgmres->modifyctx      = NULL;
733:   fgmres->modifydestroy  = NULL;
734:   fgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
735:   return(0);
736: }