Actual source code: fgmres.c

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

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

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

 12: */

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

 21: EXTERN PetscErrorCode KSPView_GMRES(KSP,PetscViewer);
 22: /*

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

 26:     This is called once, usually automatically by KSPSolveQ() or KSPSetUp(),
 27:     but can be called directly by KSPSetUp().

 29: */
 32: PetscErrorCode    KSPSetUp_FGMRES(KSP ksp)
 33: {
 34:   PetscInt       size,hh,hes,rs,cc;
 36:   PetscInt       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(PetscInt),&fgmres->mwork_alloc);
 77:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));

 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:   KSPGetVecs(ksp,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:   KSPGetVecs(ksp,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 PetscErrorCode FGMRESResidual(KSP ksp)
123: {
124:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);
125:   PetscScalar    mone = -1.0;
126:   Mat            Amat,Pmat;
127:   MatStructure   pflag;

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

133:   /* put A*x into VEC_TEMP */
134:   MatMult(Amat,ksp->vec_sol,VEC_TEMP);
135:   /* now put residual (-A*x + f) into vec_vv(0) */
136:   VecWAXPY(&mone,VEC_TEMP,ksp->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: PetscErrorCode FGMREScycle(PetscInt *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 */
171:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
172:   PetscInt       max_k = fgmres->max_k; /* max # of directions Krylov space */
173:   Mat            Amat,Pmat;
174:   MatStructure   pflag;


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

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

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

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

195:   /* FYI: AMS calls are for memory snooper */
196:   PetscObjectTakeAccess(ksp);
197:   ksp->rnorm = res_norm;
198:   PetscObjectGrantAccess(ksp);
199:   KSPLogResidualHistory(ksp,res_norm);

201:   /* check for the convergence - maybe the current guess is good enough */
202:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
203:   if (ksp->reason) {
204:     if (itcount) *itcount = 0;
205:     return(0);
206:   }

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



213: 
214:   /* MAIN ITERATION LOOP BEGINNING*/
215:   /* keep iterating until we have converged OR generated the max number
216:      of directions OR reached the max number of iterations for the method */
217:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
218:   while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
219:     KSPLogResidualHistory(ksp,res_norm);
220:     fgmres->it = (loc_it - 1);
221:     KSPMonitor(ksp,ksp->its,res_norm);

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

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

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

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

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

275: 
276:     /* Now apply rotations to new col of hessenberg (and right side of system), 
277:        calculate new rotation, and get new residual norm at the same time*/
278:     FGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
279:     if (ksp->reason) break;

281:     loc_it++;
282:     fgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
283: 
284:     PetscObjectTakeAccess(ksp);
285:     ksp->its++;
286:     ksp->rnorm = res_norm;
287:     PetscObjectGrantAccess(ksp);

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

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

301:   KSPLogResidualHistory(ksp,res_norm);

303:   /*
304:      Monitor if we know that we will not return for a restart */
305:   if (ksp->reason || ksp->its >= ksp->max_it) {
306:     KSPMonitor(ksp,ksp->its,res_norm);
307:   }

309:   if (itcount) *itcount    = loc_it;

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

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

323:   return(0);
324: }

326: /*  
327:     KSPSolve_FGMRES - This routine applies the FGMRES method.


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

333:    Output Parameter:
334: .     outits - number of iterations used

336: */

340: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
341: {
343:   PetscInt       cycle_its = 0; /* iterations done in a call to FGMREScycle */
344:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
345:   PetscTruth     diagonalscale;

348:   PCDiagonalScale(ksp->pc,&diagonalscale);
349:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);

351:   PetscObjectTakeAccess(ksp);
352:   ksp->its = 0;
353:   PetscObjectGrantAccess(ksp);

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 ksp->vec_rhs) */
359:     VecCopy(ksp->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:   while (!ksp->reason) {
366:     FGMRESResidual(ksp);
367:     if (ksp->its >= ksp->max_it) break;
368:     FGMREScycle(&cycle_its,ksp);
369:   }
370:   /* mark lack of convergence */
371:   if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

373:   return(0);
374: }

376: /*

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

380: */
383: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
384: {
385:   KSP_FGMRES     *fgmres = (KSP_FGMRES*)ksp->data;
387:   PetscInt       i;

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

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

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

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

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

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

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

431:      This is an internal routine that knows about the FGMRES internals.
432:  */
435: static PetscErrorCode BuildFgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
436: {
437:   PetscScalar    tt,zero = 0.0,one = 1.0;
439:   PetscInt       ii,k,j;
440:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)(ksp->data);

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

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

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

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

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

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

480: /*

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

485:     input parameters:

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

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

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

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

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

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

537:   /* compute new plane rotation */

539:   if (!hapend) {
540: #if defined(PETSC_USE_COMPLEX)
541:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
542: #else
543:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
544: #endif
545:     if (tt == 0.0) {
546:       ksp->reason = KSP_DIVERGED_NULL;
547:       return(0);
548:     }

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

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

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

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

578: /*

580:    FGMRESGetNewVectors - This routine allocates more work vectors, starting from 
581:                          VEC_VV(it), and more preconditioned work vectors, starting 
582:                          from PREVEC(i).

584: */
587: static PetscErrorCode FGMRESGetNewVectors(KSP ksp,PetscInt it)
588: {
589:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;
590:   PetscInt       nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
591:   PetscInt       nalloc;                      /* number to allocate */
593:   PetscInt       k;
594: 
596:   nalloc = fgmres->delta_allocate; /* number of vectors to allocate 
597:                                       in a single chunk */

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

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

608:   /* work vectors */
609:   KSPGetVecs(ksp,nalloc,&fgmres->user_work[nwork]);
610:   PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
611:   for (k=0; k < nalloc; k++) {
612:     fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
613:   }
614:   /* specify size of chunk allocated */
615:   fgmres->mwork_alloc[nwork] = nalloc;

617:   /* preconditioned vectors */
618:   KSPGetVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork]);
619:   PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
620:   for (k=0; k < nalloc; k++) {
621:     fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
622:   }

624:   /* increment the number of work vector chunks */
625:   fgmres->nwork_alloc++;
626:   return(0);
627: }

629: /* 

631:    KSPBuildSolution_FGMRES

633:      Input Parameter:
634: .     ksp - the Krylov space object
635: .     ptr-

637:    Output Parameter:
638: .     result - the solution

640:    Note: this calls BuildFgmresSoln - the same function that FGMREScycle
641:    calls directly.  

643: */
646: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
647: {
648:   KSP_FGMRES     *fgmres = (KSP_FGMRES *)ksp->data;

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


674: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
675: {
677:   PetscInt       restart,indx;
678:   PetscReal      haptol;
679:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;
680:   PetscTruth     flg;
681:   const char     *types[] = {"never","ifneeded","always"};

684:   PetscOptionsHead("KSP flexible GMRES Options");
685:     PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
686:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
687:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for declaring exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
688:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
689:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate all Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
690:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
691:     PetscOptionsLogicalGroupBegin("-ksp_gmres_classicalgramschmidt","Use classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
692:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
693:     PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Use modified Gram-Schmidt (slow but more stable)","KSPGMRESSetOrthogonalization",&flg);
694:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
695:     PetscOptionsEList("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType()",types,3,types[gmres->cgstype],&indx,&flg);
696:     if (flg) {
697:       KSPGMRESSetCGSRefinementType(ksp,(KSPGMRESCGSRefinementType)indx);
698:     }
699:     PetscOptionsName("-ksp_gmres_krylov_monitor","Graphically plot the Krylov directions","KSPSetMonitor",&flg);
700:     if (flg) {
701:       PetscViewers viewers;
702:       PetscViewersCreate(ksp->comm,&viewers);
703:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(PetscErrorCode (*)(void*))PetscViewersDestroy);
704:     }
705:     PetscOptionsLogicalGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
706:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
707:     PetscOptionsLogicalGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
708:     if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
709:   PetscOptionsTail();
710:   return(0);
711: }

713: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
714: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);

717: typedef PetscErrorCode (*FCN2)(void*);
721: PetscErrorCode KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
722: {
725:   ((KSP_FGMRES *)ksp->data)->modifypc      = fcn;
726:   ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
727:   ((KSP_FGMRES *)ksp->data)->modifyctx     = ctx;
728:   return(0);
729: }

733: EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP);
734: EXTERN PetscErrorCode KSPGMRESSetRestart_GMRES(KSP,PetscInt);
735: EXTERN PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));

740: PetscErrorCode KSPDestroy_FGMRES_Internal(KSP ksp)
741: {
742:   KSP_FGMRES     *gmres = (KSP_FGMRES*)ksp->data;
744:   PetscInt       i;

747:   /* Free the Hessenberg matrix */
748:   if (gmres->hh_origin) {
749:     PetscFree(gmres->hh_origin);
750:     gmres->hh_origin = 0;
751:   }

753:   /* Free the pointer to user variables */
754:   if (gmres->vecs) {
755:     PetscFree(gmres->vecs);
756:     gmres->vecs = 0;
757:   }
758:   if (gmres->prevecs) {
759:     PetscFree (gmres->prevecs);
760:     gmres->prevecs = 0;
761:   }

763:   /* free work vectors */
764:   for (i=0; i<gmres->nwork_alloc; i++) {
765:     VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
766:     VecDestroyVecs(gmres->prevecs_user_work[i],gmres->mwork_alloc[i]);
767:   }
768:   if (gmres->user_work)  {
769:     PetscFree(gmres->user_work);
770:     gmres->user_work = 0;
771:   }
772:   if (gmres->prevecs_user_work) {
773:     PetscFree(gmres->prevecs_user_work);
774:     gmres->prevecs_user_work = 0;
775:   }
776:   if (gmres->mwork_alloc) {
777:     PetscFree(gmres->mwork_alloc);
778:     gmres->mwork_alloc = 0;
779:   }
780:   if (gmres->nrs) {
781:     PetscFree(gmres->nrs);
782:     gmres->nrs = 0;
783:   }
784:   if (gmres->sol_temp) {
785:     VecDestroy(gmres->sol_temp);
786:     gmres->sol_temp = 0;
787:   }
788:   if (gmres->Rsvd) {
789:     PetscFree(gmres->Rsvd);
790:     gmres->Rsvd = 0;
791:   }
792:   if (gmres->Dsvd) {
793:     PetscFree(gmres->Dsvd);
794:     gmres->Dsvd = 0;
795:   }
796:   if (gmres->modifydestroy) {
797:     (*gmres->modifydestroy)(gmres->modifyctx);
798:   }

800:   gmres->vv_allocated   = 0;
801:   gmres->vecs_allocated = 0;
802:   gmres->sol_temp       = 0;
803:   gmres->nwork_alloc    = 0;
804:   return(0);
805: }

810: PetscErrorCode KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
811: {
812:   KSP_FGMRES     *gmres = (KSP_FGMRES *)ksp->data;

816:   if (max_k < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
817:   if (!ksp->setupcalled) {
818:     gmres->max_k = max_k;
819:   } else if (gmres->max_k != max_k) {
820:      gmres->max_k = max_k;
821:      ksp->setupcalled = 0;
822:      /* free the data structures, then create them again */
823:      KSPDestroy_FGMRES_Internal(ksp);
824:   }
825:   return(0);
826: }

830: EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);

833: /*MC
834:      KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.  
835:                 developed by Saad with restart


838:    Options Database Keys:
839: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
840: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
841: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 
842:                              vectors are allocated as needed)
843: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
844: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
845: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 
846:                                    stability of the classical Gram-Schmidt  orthogonalization.
847: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
848: .   -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
849: -   -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

851:    Level: beginner

853:     Notes: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations

855: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
856:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
857:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
858:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor(), KSPFGMRESSetModifyPC(),
859:            KSPFGMRESModifyPCKSP()

861: M*/

866: PetscErrorCode KSPCreate_FGMRES(KSP ksp)
867: {
868:   KSP_FGMRES     *fgmres;

872:   PetscNew(KSP_FGMRES,&fgmres);
873:   PetscMemzero(fgmres,sizeof(KSP_FGMRES));
874:   PetscLogObjectMemory(ksp,sizeof(KSP_FGMRES));
875:   ksp->data                              = (void*)fgmres;
876:   ksp->ops->buildsolution                = KSPBuildSolution_FGMRES;

878:   ksp->ops->setup                        = KSPSetUp_FGMRES;
879:   ksp->ops->solve                        = KSPSolve_FGMRES;
880:   ksp->ops->destroy                      = KSPDestroy_FGMRES;
881:   ksp->ops->view                         = KSPView_GMRES;
882:   ksp->ops->setfromoptions               = KSPSetFromOptions_FGMRES;
883:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
884:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

886:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
887:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
888:                                      KSPGMRESSetPreAllocateVectors_GMRES);
889:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
890:                                     "KSPGMRESSetOrthogonalization_GMRES",
891:                                      KSPGMRESSetOrthogonalization_GMRES);
892:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
893:                                     "KSPGMRESSetRestart_FGMRES",
894:                                      KSPGMRESSetRestart_FGMRES);
895:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
896:                                     "KSPFGMRESSetModifyPC_FGMRES",
897:                                      KSPFGMRESSetModifyPC_FGMRES);
898:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
899:                                     "KSPGMRESSetCGSRefinementType_GMRES",
900:                                      KSPGMRESSetCGSRefinementType_GMRES);


903:   fgmres->haptol              = 1.0e-30;
904:   fgmres->q_preallocate       = 0;
905:   fgmres->delta_allocate      = FGMRES_DELTA_DIRECTIONS;
906:   fgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
907:   fgmres->nrs                 = 0;
908:   fgmres->sol_temp            = 0;
909:   fgmres->max_k               = FGMRES_DEFAULT_MAXK;
910:   fgmres->Rsvd                = 0;
911:   fgmres->modifypc            = KSPFGMRESModifyPCNoChange;
912:   fgmres->modifyctx           = PETSC_NULL;
913:   fgmres->modifydestroy       = PETSC_NULL;
914:   fgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
915:   /*
916:         This is not great since it changes this without explicit request from the user
917:      but there is no left preconditioning in the FGMRES
918:   */
919:   PetscLogInfo(ksp,"Warning: Setting PC_SIDE for FGMRES to right!\n");
920:   ksp->pc_side                = PC_RIGHT;

922:   return(0);
923: }