Actual source code: fgmres.c

  1: /* $Id: fgmres.c,v 1.29 2001/08/07 21:30:49 bsmith Exp $ */

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

  7:     Preconditioning:  It the preconditioner is constant then this fgmres
  8:     code is equivalent to RIGHT-PRECONDITIONED GMRES.

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

 14: */

 16:  #include src/sles/ksp/impls/fgmres/fgmresp.h
 17: #define FGMRES_DELTA_DIRECTIONS 10
 18: #define FGMRES_DEFAULT_MAXK     30
 19: static int    FGMRESGetNewVectors(KSP,int);
 20: static int    FGMRESUpdateHessenberg(KSP,int,PetscTruth,PetscReal *);
 21: static int    BuildFgmresSoln(PetscScalar*,Vec,Vec,KSP,int);

 23: /*

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

 27:     This is called once, usually automatically by SLESSolve() or SLESSetUp(),
 28:     but can be called directly by KSPSetUp().

 30: */
 33: int    KSPSetUp_FGMRES(KSP ksp)
 34: {
 35:   unsigned  int size,hh,hes,rs,cc;
 36:   int           ierr,max_k,k;
 37:   KSP_FGMRES    *fgmres = (KSP_FGMRES *)ksp->data;

 40:   if (ksp->pc_side == PC_SYMMETRIC) {
 41:     SETERRQ(2,"no symmetric preconditioning for KSPFGMRES");
 42:   } else if (ksp->pc_side == PC_LEFT) {
 43:     SETERRQ(2,"no left preconditioning for KSPFGMRES");
 44:   }
 45:   max_k         = fgmres->max_k;
 46:   hh            = (max_k + 2) * (max_k + 1);
 47:   hes           = (max_k + 1) * (max_k + 1);
 48:   rs            = (max_k + 2);
 49:   cc            = (max_k + 1);  /* SS and CC are the same size */
 50:   size          = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);

 52:   /* Allocate space and set pointers to beginning */
 53:   PetscMalloc(size,&fgmres->hh_origin);
 54:   PetscMemzero(fgmres->hh_origin,size);
 55:   PetscLogObjectMemory(ksp,size);                      /* HH - modified (by plane 
 56:                                                      rotations) hessenburg */
 57:   fgmres->hes_origin = fgmres->hh_origin + hh;     /* HES - unmodified hessenburg */
 58:   fgmres->rs_origin  = fgmres->hes_origin + hes;   /* RS - the right-hand-side of the 
 59:                                                       Hessenberg system */
 60:   fgmres->cc_origin  = fgmres->rs_origin + rs;     /* CC - cosines for rotations */
 61:   fgmres->ss_origin  = fgmres->cc_origin + cc;     /* SS - sines for rotations */

 63:   if (ksp->calc_sings) {
 64:     /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
 65:     size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
 66:     PetscMalloc(size,&fgmres->Rsvd);
 67:     PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&fgmres->Dsvd);
 68:     PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
 69:   }

 71:   /* Allocate array to hold pointers to user vectors.  Note that we need
 72:    4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
 73:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&fgmres->vecs);
 74:   fgmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
 75:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&fgmres->user_work);
 76:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(int),&fgmres->mwork_alloc);
 77:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)+sizeof(int)));

 79:   /* New for FGMRES - Allocate array to hold pointers to preconditioned 
 80:      vectors - same sizes as user vectors above */
 81:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&fgmres->prevecs);
 82:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&fgmres->prevecs_user_work);
 83:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)));


 86:   /* if q_preallocate = 0 then only allocate one "chunck" of space (for 
 87:      5 vectors) - additional will then be allocated from FGMREScycle() 
 88:      as needed.  Otherwise, allocate all of the space that could be needed */
 89:   if (fgmres->q_preallocate) {
 90:     fgmres->vv_allocated   = VEC_OFFSET + 2 + max_k;
 91:   } else {
 92:     fgmres->vv_allocated    = 5;
 93:   }

 95:   /* space for work vectors */
 96:   VecDuplicateVecs(VEC_RHS,fgmres->vv_allocated,&fgmres->user_work[0]);
 97:   PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->user_work[0]);
 98:   for (k=0; k < fgmres->vv_allocated; k++) {
 99:     fgmres->vecs[k] = fgmres->user_work[0][k];
100:   }

102:   /* space for preconditioned vectors */
103:   VecDuplicateVecs(VEC_RHS,fgmres->vv_allocated,&fgmres->prevecs_user_work[0]);
104:   PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->prevecs_user_work[0]);
105:   for (k=0; k < fgmres->vv_allocated; k++) {
106:     fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
107:   }

109:   /* specify how many work vectors have been allocated in this 
110:      chunck" (the first one) */
111:   fgmres->mwork_alloc[0] = fgmres->vv_allocated;
112:   fgmres->nwork_alloc    = 1;

114:   return(0);
115: }

117: /* 
118:     FGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED) 
119: */
122: static int FGMRESResidual(KSP ksp)
123: {
124:   KSP_FGMRES   *fgmres = (KSP_FGMRES *)(ksp->data);
125:   PetscScalar  mone = -1.0;
126:   Mat          Amat,Pmat;
127:   MatStructure pflag;
128:   int          ierr;

131:   PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);

133:   /* put A*x into VEC_TEMP */
134:   MatMult(Amat,VEC_SOLN,VEC_TEMP);
135:   /* now put residual (-A*x + f) into vec_vv(0) */
136:   VecWAXPY(&mone,VEC_TEMP,VEC_RHS,VEC_VV(0));
137:   return(0);
138: }

140: /*

142:     FGMRESCycle - Run fgmres, possibly with restart.  Return residual 
143:                   history if requested.

145:     input parameters:
146: .         fgmres  - structure containing parameters and work areas

148:     output parameters:
149: .        itcount - number of iterations used.  If null, ignored.
150: .        converged - 0 if not converged

152:                   
153:     Notes:
154:     On entry, the value in vector VEC_VV(0) should be 
155:     the initial residual.


158:  */
161: int FGMREScycle(int *itcount,KSP ksp)
162: {

164:   KSP_FGMRES   *fgmres = (KSP_FGMRES *)(ksp->data);
165:   PetscReal    res_norm;
166:   PetscReal    hapbnd,tt;
167:   PetscScalar  zero = 0.0;
168:   PetscScalar  tmp;
169:   PetscTruth   hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
170:   int          ierr;
171:   int          loc_it;                /* local count of # of dir. in Krylov space */
172:   int          max_k = fgmres->max_k; /* max # of directions Krylov space */
173:   int          max_it = ksp->max_it;  /* max # of overall iterations for the method */
174:   Mat          Amat,Pmat;
175:   MatStructure pflag;


179:   /* Number of pseudo iterations since last restart is the number 
180:      of prestart directions */
181:   loc_it = 0;

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

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

190:   /* check for the convergence - maybe the current guess is good enough */
191:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
192:   if (ksp->reason) {
193:     if (itcount) *itcount = 0;
194:     return(0);
195:   }

197:   /* scale VEC_VV (the initial residual) */
198:   tmp = 1.0/res_norm; VecScale(&tmp,VEC_VV(0));

200:   /* FYI: AMS calls are for memory snooper */
201:   PetscObjectTakeAccess(ksp);
202:   ksp->rnorm = res_norm;
203:   PetscObjectGrantAccess(ksp);


206:   /* note: (fgmres->it) is always set one less than (loc_it) It is used in 
207:      KSPBUILDSolution_FGMRES, where it is passed to BuildFGmresSoln.  
208:      Note that when BuildFGmresSoln is called from this function, 
209:      (loc_it -1) is passed, so the two are equivalent */
210:   fgmres->it = (loc_it - 1);
211: 
212:   /* MAIN ITERATION LOOP BEGINNING*/
213:   /* keep iterating until we have converged OR generated the max number
214:      of directions OR reached the max number of iterations for the method */
215:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
216:   while (!ksp->reason && loc_it < max_k && ksp->its < max_it) {
217:     KSPLogResidualHistory(ksp,res_norm);
218:     fgmres->it = (loc_it - 1);
219:     KSPMonitor(ksp,ksp->its,res_norm);

221:     /* see if more space is needed for work vectors */
222:     if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
223:       FGMRESGetNewVectors(ksp,loc_it+1);
224:       /* (loc_it+1) is passed in as number of the first vector that should
225:          be allocated */
226:     }

228:     /* CHANGE THE PRECONDITIONER? */
229:     /* ModifyPC is the callback function that can be used to
230:        change the PC or its attributes before its applied */
231:     (*fgmres->modifypc)(ksp,ksp->its,loc_it,res_norm,fgmres->modifyctx);
232: 
233: 
234:     /* apply PRECONDITIONER to direction vector and store with 
235:        preconditioned vectors in prevec */
236:     PCApply(ksp->B,VEC_VV(loc_it),PREVEC(loc_it),PC_RIGHT);
237: 
238:     PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
239:     /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
240:     MatMult(Amat,PREVEC(loc_it),VEC_VV(1+loc_it));

242: 
243:     /* update hessenberg matrix and do Gram-Schmidt - new direction is in
244:        VEC_VV(1+loc_it)*/
245:     (*fgmres->orthog)(ksp,loc_it);

247:     /* new entry in hessenburg is the 2-norm of our new direction */
248:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
249:     *HH(loc_it+1,loc_it)   = tt;
250:     *HES(loc_it+1,loc_it)  = tt;

252:     /* Happy Breakdown Check */
253:     hapbnd  = PetscAbsScalar((tt) / *RS(loc_it));
254:     /* RS(loc_it) contains the res_norm from the last iteration  */
255:     hapbnd = PetscMin(fgmres->haptol,hapbnd);
256:     if (tt > hapbnd) {
257:         tmp = 1.0/tt;
258:         /* scale new direction by its norm */
259:         VecScale(&tmp,VEC_VV(loc_it+1));
260:     } else {
261:         /* This happens when the solution is exactly reached. */
262:         /* So there is no new direction... */
263:           VecSet(&zero,VEC_TEMP); /* set VEC_TEMP to 0 */
264:           hapend = PETSC_TRUE;
265:     }
266:     /* note that for FGMRES we could get HES(loc_it+1, loc_it)  = 0 and the
267:        current solution would not be exact if HES was singular.  Note that 
268:        HH non-singular implies that HES is no singular, and HES is guaranteed
269:        to be nonsingular when PREVECS are linearly independent and A is 
270:        nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity 
271:        of HES). So we should really add a check to verify that HES is nonsingular.*/

273: 
274:     /* Now apply rotations to new col of hessenberg (and right side of system), 
275:        calculate new rotation, and get new residual norm at the same time*/
276:     FGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
277:     loc_it++;
278:     fgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
279: 
280:     PetscObjectTakeAccess(ksp);
281:     ksp->its++;
282:     ksp->rnorm = res_norm;
283:     PetscObjectGrantAccess(ksp);

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

287:     /* Catch error in happy breakdown and signal convergence and break from loop */
288:     if (hapend) {
289:       if (!ksp->reason) {
290:         SETERRQ(0,"You reached the happy break down,but convergence was not indicated.");
291:       }
292:       break;
293:     }
294:   }
295:   /* END OF ITERATION LOOP */

297:   KSPLogResidualHistory(ksp,res_norm);

299:   /*
300:      Monitor if we know that we will not return for a restart */
301:   if (ksp->reason || ksp->its >= max_it) {
302:     KSPMonitor(ksp,ksp->its,res_norm);
303:   }

305:   if (itcount) *itcount    = loc_it;

307:   /*
308:     Down here we have to solve for the "best" coefficients of the Krylov
309:     columns, add the solution values together, and possibly unwind the
310:     preconditioning from the solution
311:    */
312: 
313:   /* Form the solution (or the solution so far) */
314:   /* Note: must pass in (loc_it-1) for iteration count so that BuildFgmresSoln
315:      properly navigates */

317:   BuildFgmresSoln(RS(0),VEC_SOLN,VEC_SOLN,ksp,loc_it-1);

319:   return(0);
320: }

322: /*  
323:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

329:    Output Parameter:
330: .     outits - number of iterations used

332: */

336: int KSPSolve_FGMRES(KSP ksp,int *outits)
337: {
338:   int        ierr;
339:   int        cycle_its; /* iterations done in a call to FGMREScycle */
340:   int        itcount;   /* running total of iterations, incl. those in restarts */
341:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
342:   PetscTruth diagonalscale;

345:   PCDiagonalScale(ksp->B,&diagonalscale);
346:   if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);

348:   PetscObjectTakeAccess(ksp);
349:   ksp->its = 0;
350:   PetscObjectGrantAccess(ksp);

352:   /* initialize */
353:   itcount  = 0;

355:   /* Compute the initial (NOT preconditioned) residual */
356:   if (!ksp->guess_zero) {
357:     FGMRESResidual(ksp);
358:   } else { /* guess is 0 so residual is F (which is in VEC_RHS) */
359:     VecCopy(VEC_RHS,VEC_VV(0));
360:   }
361:   /* now the residual is in VEC_VV(0) - which is what 
362:      FGMREScycle expects... */
363: 
364:   FGMREScycle(&cycle_its,ksp);
365:   itcount += cycle_its;
366:   while (!ksp->reason) {
367:     FGMRESResidual(ksp);
368:     if (itcount >= ksp->max_it) break;
369:     FGMREScycle(&cycle_its,ksp);
370:     itcount += cycle_its;
371:   }
372:   /* mark lack of convergence */
373:   if (itcount >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

375:   *outits = itcount;
376:   return(0);
377: }

379: /*

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

383: */
386: int KSPDestroy_FGMRES(KSP ksp)
387: {
388:   KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
389:   int       i,ierr;

392:   /* Free the Hessenberg matrices */
393:   if (fgmres->hh_origin) {PetscFree(fgmres->hh_origin);}

395:   /* Free pointers to user variables */
396:   if (fgmres->vecs) {PetscFree(fgmres->vecs);}
397:   if (fgmres->prevecs) {PetscFree (fgmres->prevecs);}

399:   /* free work vectors */
400:   for (i=0; i < fgmres->nwork_alloc; i++) {
401:     VecDestroyVecs(fgmres->user_work[i],fgmres->mwork_alloc[i]);
402:   }
403:   if (fgmres->user_work)  {PetscFree(fgmres->user_work);}

405:   for (i=0; i < fgmres->nwork_alloc; i++) {
406:     VecDestroyVecs(fgmres->prevecs_user_work[i],fgmres->mwork_alloc[i]);
407:   }
408:   if (fgmres->prevecs_user_work) {PetscFree(fgmres->prevecs_user_work);}

410:   if (fgmres->mwork_alloc) {PetscFree(fgmres->mwork_alloc);}
411:   if (fgmres->nrs) {PetscFree(fgmres->nrs);}
412:   if (fgmres->sol_temp) {VecDestroy(fgmres->sol_temp);}
413:   if (fgmres->Rsvd) {PetscFree(fgmres->Rsvd);}
414:   if (fgmres->Dsvd) {PetscFree(fgmres->Dsvd);}
415:   if (fgmres->modifydestroy) {
416:     (*fgmres->modifydestroy)(fgmres->modifyctx);
417:   }
418:   PetscFree(fgmres);
419:   return(0);
420: }

422: /*
423:     BuildFgmresSoln - create the solution from the starting vector and the
424:                       current iterates.

426:     Input parameters:
427:         nrs - work area of size it + 1.
428:         vguess  - index of initial guess
429:         vdest - index of result.  Note that vguess may == vdest (replace
430:                 guess with the solution).
431:         it - HH upper triangular part is a block of size (it+1) x (it+1)  

433:      This is an internal routine that knows about the FGMRES internals.
434:  */
437: static int BuildFgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,int it)
438: {
439:   PetscScalar  tt,zero = 0.0,one = 1.0;
440:   int          ierr,ii,k,j;
441:   KSP_FGMRES   *fgmres = (KSP_FGMRES *)(ksp->data);

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

446:   /* If it is < 0, no fgmres steps have been performed */
447:   if (it < 0) {
448:     if (vdest != vguess) {
449:       VecCopy(vguess,vdest);
450:     }
451:     return(0);
452:   }

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

456:   /* solve the upper triangular system - RS is the right side and HH is 
457:      the upper triangular matrix  - put soln in nrs */
458:   nrs[it] = *RS(it) / *HH(it,it);
459:   for (ii=1; ii<=it; ii++) {
460:     k   = it - ii;
461:     tt  = *RS(k);
462:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
463:     nrs[k]   = tt / *HH(k,k);
464:   }

466:   /* Accumulate the correction to the soln of the preconditioned prob. in 
467:      VEC_TEMP - note that we use the preconditioned vectors  */
468:   VecSet(&zero,VEC_TEMP); /* set VEC_TEMP components to 0 */
469:   VecMAXPY(it+1,nrs,VEC_TEMP,&PREVEC(0));

471:   /* put updated solution into vdest.*/
472:   if (vdest != vguess) {
473:     VecCopy(VEC_TEMP,vdest);
474:     VecAXPY(&one,vguess,vdest);
475:   } else  {/* replace guess with solution */
476:     VecAXPY(&one,VEC_TEMP,vdest);
477:   }
478:   return(0);
479: }

481: /*

483:     FGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
484:                             Return new residual.

486:     input parameters:

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

493:     output parameters:
494: .        res - the new residual
495:         
496:  */
499: static int FGMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
500: {
501:   PetscScalar   *hh,*cc,*ss,tt;
502:   int           j;
503:   KSP_FGMRES    *fgmres = (KSP_FGMRES *)(ksp->data);

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

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

516:   for (j=1; j<=it; j++) {
517:     tt  = *hh;
518: #if defined(PETSC_USE_COMPLEX)
519:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
520: #else
521:     *hh = *cc * tt + *ss * *(hh+1);
522: #endif
523:     hh++;
524:     *hh = *cc++ * *hh - (*ss++ * tt);
525:     /* hh, cc, and ss have all been incremented one by end of loop */
526:   }

528:   /*
529:     compute the new plane rotation, and apply it to:
530:      1) the right-hand-side of the Hessenberg system (RS)
531:         note: it affects RS(it) and RS(it+1)
532:      2) the new column of the Hessenberg matrix
533:         note: it affects HH(it,it) which is currently pointed to 
534:         by hh and HH(it+1, it) (*(hh+1))  
535:     thus obtaining the updated value of the residual...
536:   */

538:   /* compute new plane rotation */

540:   if (!hapend) {
541: #if defined(PETSC_USE_COMPLEX)
542:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
543: #else
544:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
545: #endif
546:     if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
547:     *cc       = *hh / tt;   /* new cosine value */
548:     *ss       = *(hh+1) / tt;  /* new sine value */

550:     /* apply to 1) and 2) */
551:     *RS(it+1) = - (*ss * *RS(it));
552: #if defined(PETSC_USE_COMPLEX)
553:     *RS(it)   = PetscConj(*cc) * *RS(it);
554:     *hh       = PetscConj(*cc) * *hh + *ss * *(hh+1);
555: #else
556:     *RS(it)   = *cc * *RS(it);
557:     *hh       = *cc * *hh + *ss * *(hh+1);
558: #endif

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

563:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
564:             another rotation matrix (so RH doesn't change).  The new residual is 
565:             always the new sine term times the residual from last time (RS(it)), 
566:             but now the new sine rotation would be zero...so the residual should
567:             be zero...so we will multiply "zero" by the last residual.  This might
568:             not be exactly what we want to do here -could just return "zero". */
569: 
570:     *res = 0.0;
571:   }
572:   return(0);
573: }

575: /*

577:    FGMRESGetNewVectors - This routine allocates more work vectors, starting from 
578:                          VEC_VV(it), and more preconditioned work vectors, starting 
579:                          from PREVEC(i).

581: */
584: static int FGMRESGetNewVectors(KSP ksp,int it)
585: {
586:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
587:   int        nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
588:   int        nalloc;                      /* number to allocate */
589:   int        k,ierr;
590: 
592:   nalloc = fgmres->delta_allocate; /* number of vectors to allocate 
593:                                       in a single chunk */

595:   /* Adjust the number to allocate to make sure that we don't exceed the
596:      number of available slots (fgmres->vecs_allocated)*/
597:   if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated){
598:     nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
599:   }
600:   if (!nalloc) return(0);

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

604:   /* work vectors */
605:   VecDuplicateVecs(VEC_RHS,nalloc,&fgmres->user_work[nwork]);
606:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
607:   for (k=0; k < nalloc; k++) {
608:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
609:   }
610:   /* specify size of chunk allocated */
611:   fgmres->mwork_alloc[nwork] = nalloc;

613:   /* preconditioned vectors */
614:   VecDuplicateVecs(VEC_RHS,nalloc,&fgmres->prevecs_user_work[nwork]);
615:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
616:   for (k=0; k < nalloc; k++) {
617:     fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
618:   }

620:   /* increment the number of work vector chunks */
621:   fgmres->nwork_alloc++;
622:   return(0);
623: }

625: /* 

627:    KSPBuildSolution_FGMRES

629:      Input Parameter:
630: .     ksp - the Krylov space object
631: .     ptr-

633:    Output Parameter:
634: .     result - the solution

636:    Note: this calls BuildFgmresSoln - the same function that FGMREScycle
637:    calls directly.  

639: */
642: int KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
643: {
644:   KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
645:   int        ierr;

648:   if (!ptr) {
649:     if (!fgmres->sol_temp) {
650:       VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
651:       PetscLogObjectParent(ksp,fgmres->sol_temp);
652:     }
653:     ptr = fgmres->sol_temp;
654:   }
655:   if (!fgmres->nrs) {
656:     /* allocate the work area */
657:     PetscMalloc(fgmres->max_k*sizeof(PetscScalar),&fgmres->nrs);
658:     PetscLogObjectMemory(ksp,fgmres->max_k*sizeof(PetscScalar));
659:   }
660: 
661:   BuildFgmresSoln(fgmres->nrs,VEC_SOLN,ptr,ksp,fgmres->it);
662:   *result = ptr;
663: 
664:   return(0);
665: }

667: /*

669:    KSPView_FGMRES -Prints information about the current Krylov method 
670:                   being used.

672:  */
675: int KSPView_FGMRES(KSP ksp,PetscViewer viewer)
676: {
677:   KSP_FGMRES   *fgmres = (KSP_FGMRES *)ksp->data;
678:   char         *cstr;
679:   int          ierr;
680:   PetscTruth   isascii;

683:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
684:   if (isascii) {
685:     if (fgmres->orthog == KSPGMRESUnmodifiedGramSchmidtOrthogonalization) {
686:       cstr = "Unmodified Gram-Schmidt Orthogonalization";
687:     } else if (fgmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
688:       cstr = "Modified Gram-Schmidt Orthogonalization";
689:     } else if (fgmres->orthog == KSPGMRESIROrthogonalization) {
690:       cstr = "Unmodified Gram-Schmidt + 1 step Iterative Refinement Orthogonalization";
691:     } else {
692:       cstr = "unknown orthogonalization";
693:     }
694:     PetscViewerASCIIPrintf(viewer,"  FGMRES: restart=%d, using %s\n",fgmres->max_k,cstr);
695:   } else {
696:     SETERRQ(1,"Viewer type not supported for this object");
697:   }
698:   return(0);
699: }

703: int KSPSetFromOptions_FGMRES(KSP ksp)
704: {
705:   int         ierr,restart;
706:   PetscReal   haptol;
707:   KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;
708:   PetscTruth  flg;

711:   PetscOptionsHead("KSP flexible GMRES Options");
712:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
713:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
714:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for declaring exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
715:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
716:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate all Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
717:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
718:     PetscOptionsLogicalGroupBegin("-ksp_gmres_unmodifiedgramschmidt","Use classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
719:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESUnmodifiedGramSchmidtOrthogonalization);}
720:     PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Use modified Gram-Schmidt (slow but more stable)","KSPGMRESSetOrthogonalization",&flg);
721:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
722:     PetscOptionsLogicalGroupEnd("-ksp_gmres_irorthog","Use classical Gram-Schmidt with iterative refinement","KSPGMRESSetOrthogonalization",&flg);
723:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESIROrthogonalization);}
724:     PetscOptionsName("-ksp_gmres_krylov_monitor","Graphically plot the Krylov directions","KSPSetMonitor",&flg);
725:     if (flg) {
726:       PetscViewers viewers;
727:       PetscViewersCreate(ksp->comm,&viewers);
728:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
729:     }
730:     PetscOptionsLogicalGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
731:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
732:     PetscOptionsLogicalGroupEnd("-ksp_fgmres_modifypcsles","vary the SLES based preconditioner","KSPFGMRESSetModifyPC",&flg);
733:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCSLES,0,0);}
734:   PetscOptionsTail();
735:   return(0);
736: }

738: EXTERN int KSPComputeExtremeSingularValues_FGMRES(KSP,PetscReal *,PetscReal *);
739: EXTERN int KSPComputeEigenvalues_FGMRES(KSP,int,PetscReal *,PetscReal *,int *);

741: EXTERN_C_BEGIN
744: int KSPFGMRESSetModifyPC_FGMRES(KSP ksp,int (*fcn)(KSP,int,int,PetscReal,void*),void *ctx,int (*d)(void*))
745: {
748:   ((KSP_FGMRES *)ksp->data)->modifypc      = fcn;
749:   ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
750:   ((KSP_FGMRES *)ksp->data)->modifyctx     = ctx;
751:   return(0);
752: }
753: EXTERN_C_END

755: EXTERN_C_BEGIN
756: EXTERN int KSPGMRESSetPreAllocateVectors_GMRES(KSP);
757: EXTERN int KSPGMRESSetRestart_GMRES(KSP,int);
758: EXTERN int KSPGMRESSetOrthogonalization_GMRES(KSP,int (*)(KSP,int));
759: EXTERN_C_END

763: int KSPDestroy_FGMRES_Internal(KSP ksp)
764: {
765:   KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;
766:   int       i,ierr;

769:   /* Free the Hessenberg matrix */
770:   if (gmres->hh_origin) {PetscFree(gmres->hh_origin);}

772:   /* Free the pointer to user variables */
773:   if (gmres->vecs) {PetscFree(gmres->vecs);}

775:   /* free work vectors */
776:   for (i=0; i<gmres->nwork_alloc; i++) {
777:     VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
778:   }
779:   if (gmres->user_work)  {PetscFree(gmres->user_work);}
780:   if (gmres->mwork_alloc) {PetscFree(gmres->mwork_alloc);}
781:   if (gmres->nrs) {PetscFree(gmres->nrs);}
782:   if (gmres->sol_temp) {VecDestroy(gmres->sol_temp);}
783:   if (gmres->Rsvd) {PetscFree(gmres->Rsvd);}
784:   if (gmres->Dsvd) {PetscFree(gmres->Dsvd);}

786:   return(0);
787: }

789: EXTERN_C_BEGIN
792: int KSPGMRESSetRestart_FGMRES(KSP ksp,int max_k)
793: {
794:   KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;
795:   int        ierr;

798:   if (max_k < 1) SETERRQ(1,"Restart must be positive");
799:   if (!ksp->setupcalled) {
800:     gmres->max_k = max_k;
801:   } else if (gmres->max_k != max_k) {
802:      gmres->max_k = max_k;
803:      ksp->setupcalled = 0;
804:      /* free the data structures, then create them again */
805:      KSPDestroy_FGMRES_Internal(ksp);
806:   }
807:   return(0);
808: }
809: EXTERN_C_END

811: EXTERN_C_BEGIN
814: int KSPCreate_FGMRES(KSP ksp)
815: {
816:   KSP_FGMRES *fgmres;
817:   int        ierr;

820:   PetscNew(KSP_FGMRES,&fgmres);
821:   PetscMemzero(fgmres,sizeof(KSP_FGMRES));
822:   PetscLogObjectMemory(ksp,sizeof(KSP_FGMRES));
823:   ksp->data                              = (void*)fgmres;
824:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;

826:   ksp->ops->setup                        = KSPSetUp_FGMRES;
827:   ksp->ops->solve                        = KSPSolve_FGMRES;
828:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
829:   ksp->ops->view                         = KSPView_FGMRES;
830:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
831:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FGMRES;
832:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_FGMRES;

834:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
835:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
836:                                      KSPGMRESSetPreAllocateVectors_GMRES);
837:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
838:                                     "KSPGMRESSetOrthogonalization_GMRES",
839:                                      KSPGMRESSetOrthogonalization_GMRES);
840:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
841:                                     "KSPGMRESSetRestart_FGMRES",
842:                                      KSPGMRESSetRestart_FGMRES);
843:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
844:                                     "KSPFGMRESSetModifyPC_FGMRES",
845:                                      KSPFGMRESSetModifyPC_FGMRES);


848:   fgmres->haptol              = 1.0e-30;
849:   fgmres->q_preallocate       = 0;
850:   fgmres->delta_allocate      = FGMRES_DELTA_DIRECTIONS;
851:   fgmres->orthog              = KSPGMRESUnmodifiedGramSchmidtOrthogonalization;
852:   fgmres->nrs                 = 0;
853:   fgmres->sol_temp            = 0;
854:   fgmres->max_k               = FGMRES_DEFAULT_MAXK;
855:   fgmres->Rsvd                = 0;
856:   fgmres->modifypc            = KSPFGMRESModifyPCNoChange;
857:   fgmres->modifyctx           = PETSC_NULL;
858:   fgmres->modifydestroy       = PETSC_NULL;

860:   /*
861:         This is not great since it changes this without explicit request from the user
862:      but there is no left preconditioning in the FGMRES
863:   */
864:   PetscLogInfo(ksp,"Warning: Setting PC_SIDE for FGMRES to right!\n");
865:   ksp->pc_side                = PC_RIGHT;

867:   return(0);
868: }
869: EXTERN_C_END