Actual source code: fgmres.c

petsc-3.4.2 2013-07-02
  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);

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

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

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

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

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


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

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


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

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

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

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


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

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

213:     PetscObjectAMSTakeAccess((PetscObject)ksp);
214:     ksp->its++;
215:     ksp->rnorm = res_norm;
216:     PetscObjectAMSGrantAccess((PetscObject)ksp);

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

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

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

240:   if (itcount) *itcount = loc_it;

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

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

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

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


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

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

266: */

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

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

281:   PetscObjectAMSTakeAccess((PetscObject)ksp);
282:   ksp->its = 0;
283:   PetscObjectAMSGrantAccess((PetscObject)ksp);

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

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

305: extern PetscErrorCode KSPReset_FGMRES(KSP);
306: /*

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

310: */
313: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
314: {

318:   KSPReset_FGMRES(ksp);
319:   PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",NULL);
320:   KSPDestroy_GMRES(ksp);
321:   return(0);
322: }

324: /*
325:     KSPFGMRESBuildSoln - create the solution from the starting vector and the
326:                       current iterates.

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

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

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

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

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

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

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

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

386: /*

388:     KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
389:                             Return new residual.

391:     input parameters:

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

398:     output parameters:
399: .        res - the new residual

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

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

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

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

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

439:   /* compute new plane rotation */

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

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

451:     /* apply to 1) and 2) */
452:     *RS(it+1) = -(*ss * *RS(it));
453:     *RS(it)   = PetscConj(*cc) * *RS(it);
454:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);

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

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

466:     *res = 0.0;
467:   }
468:   return(0);
469: }

471: /*

473:    KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
474:                          VEC_VV(it), and more preconditioned work vectors, starting
475:                          from PREVEC(i).

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

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

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

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

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

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

517:   /* increment the number of work vector chunks */
518:   fgmres->nwork_alloc++;
519:   return(0);
520: }

522: /*

524:    KSPBuildSolution_FGMRES

526:      Input Parameter:
527: .     ksp - the Krylov space object
528: .     ptr-

530:    Output Parameter:
531: .     result - the solution

533:    Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
534:    calls directly.

536: */
539: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
540: {
541:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;

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

558:   KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
559:   if (result) *result = ptr;
560:   return(0);
561: }

565: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
566: {
568:   PetscBool      flg;

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

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

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


599: PetscErrorCode KSPReset_FGMRES(KSP ksp)
600: {
601:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
603:   PetscInt       i;

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

620: PetscErrorCode  KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
621: {
622:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;

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

640: PetscErrorCode  KSPGMRESGetRestart_FGMRES(KSP ksp,PetscInt *max_k)
641: {
642:   KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;

645:   *max_k = gmres->max_k;
646:   return(0);
647: }

649: /*MC
650:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
651:                 developed by Saad with restart


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

667:    Level: beginner

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

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

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

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

683: M*/

687: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
688: {
689:   KSP_FGMRES     *fgmres;

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

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

706:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
707:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,0);

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


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