Actual source code: lgmres.c

  1: #define PETSCKSP_DLL

 3:  #include ../src/ksp/ksp/impls/gmres/lgmres/lgmresimpl.h

  5: #define LGMRES_DELTA_DIRECTIONS 10
  6: #define LGMRES_DEFAULT_MAXK     30
  7: #define LGMRES_DEFAULT_AUGDIM   2 /*default number of augmentation vectors */ 
  8: static PetscErrorCode    LGMRESGetNewVectors(KSP,PetscInt);
  9: static PetscErrorCode    LGMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal *);
 10: static PetscErrorCode    BuildLgmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 14: PetscErrorCode  KSPLGMRESSetAugDim(KSP ksp, PetscInt dim)
 15: {

 19:   PetscTryMethod((ksp),"KSPLGMRESSetAugDim_C",(KSP,PetscInt),(ksp,dim));
 20:   return(0);
 21: }

 25: PetscErrorCode  KSPLGMRESSetConstant(KSP ksp)
 26: {

 30:   PetscTryMethod((ksp),"KSPLGMRESSetConstant_C",(KSP),(ksp));
 31:   return(0);
 32: }

 35: /*
 36:     KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.

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

 41: */
 44: PetscErrorCode    KSPSetUp_LGMRES(KSP ksp)
 45: {
 47:   PetscInt       max_k,k, aug_dim;
 48:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;

 51:   if (ksp->pc_side == PC_SYMMETRIC) {
 52:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPLGMRES");
 53:   }
 54:   max_k         = lgmres->max_k;
 55:   aug_dim       = lgmres->aug_dim;
 56:   KSPSetUp_GMRES(ksp);

 58:   /* need array of pointers to augvecs*/
 59:   PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs);
 60:   lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
 61:   PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs_user_work);
 62:   PetscMalloc(aug_dim*sizeof(PetscInt),&lgmres->aug_order);
 63:   PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));

 65:   /*  for now we will preallocate the augvecs - because aug_dim << restart
 66:      ... also keep in mind that we need to keep augvecs from cycle to cycle*/
 67:   lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
 68:   lgmres->augwork_alloc =  2* aug_dim + AUG_OFFSET;
 69:   KSPGetVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,PETSC_NULL);
 70:   PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
 71:   for (k=0; k<lgmres->aug_vv_allocated; k++) {
 72:     lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
 73:   }
 74:   return(0);
 75: }


 78: /*

 80:     LGMRESCycle - Run lgmres, possibly with restart.  Return residual 
 81:                   history if requested.

 83:     input parameters:
 84: .         lgmres  - structure containing parameters and work areas

 86:     output parameters:
 87: .        nres    - residuals (from preconditioned system) at each step.
 88:                   If restarting, consider passing nres+it.  If null, 
 89:                   ignored
 90: .        itcount - number of iterations used.   nres[0] to nres[itcount]
 91:                   are defined.  If null, ignored.  If null, ignored.
 92: .        converged - 0 if not converged

 94:                   
 95:     Notes:
 96:     On entry, the value in vector VEC_VV(0) should be 
 97:     the initial residual.


100:  */
103: PetscErrorCode LGMREScycle(PetscInt *itcount,KSP ksp)
104: {
105:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)(ksp->data);
106:   PetscReal      res_norm, res;
107:   PetscReal      hapbnd, tt;
108:   PetscScalar    tmp;
109:   PetscTruth     hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
111:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
112:   PetscInt       max_k = lgmres->max_k; /* max approx space size */
113:   PetscInt       max_it = ksp->max_it;  /* max # of overall iterations for the method */
114:   /* LGMRES_MOD - new variables*/
115:   PetscInt       aug_dim = lgmres->aug_dim;
116:   PetscInt       spot = 0;
117:   PetscInt       order = 0;
118:   PetscInt       it_arnoldi;             /* number of arnoldi steps to take */
119:   PetscInt       it_total;               /* total number of its to take (=approx space size)*/
120:   PetscInt       ii, jj;
121:   PetscReal      tmp_norm;
122:   PetscScalar    inv_tmp_norm;
123:   PetscScalar    *avec;

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

130:   /* LGMRES_MOD: determine number of arnoldi steps to take */
131:   /* if approx_constant then we keep the space the same size even if 
132:      we don't have the full number of aug vectors yet*/
133:   if (lgmres->approx_constant) {
134:      it_arnoldi = max_k - lgmres->aug_ct;
135:   } else {
136:       it_arnoldi = max_k - aug_dim;
137:   }

139:   it_total =  it_arnoldi + lgmres->aug_ct;

141:   /* initial residual is in VEC_VV(0)  - compute its norm*/
142:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
143:   res    = res_norm;
144: 
145:   /* first entry in right-hand-side of hessenberg system is just 
146:      the initial residual norm */
147:   *GRS(0) = res_norm;

149:  /* check for the convergence */
150:   if (!res) {
151:      if (itcount) *itcount = 0;
152:      ksp->reason = KSP_CONVERGED_ATOL;
153:      PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
154:      return(0);
155:   }

157:   /* scale VEC_VV (the initial residual) */
158:   tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);

160:   ksp->rnorm = res;


163:   /* note: (lgmres->it) is always set one less than (loc_it) It is used in 
164:      KSPBUILDSolution_LGMRES, where it is passed to BuildLgmresSoln.  
165:      Note that when BuildLgmresSoln is called from this function, 
166:      (loc_it -1) is passed, so the two are equivalent */
167:   lgmres->it = (loc_it - 1);

169: 
170:   /* MAIN ITERATION LOOP BEGINNING*/


173:   /* keep iterating until we have converged OR generated the max number
174:      of directions OR reached the max number of iterations for the method */
175:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
176: 
177:   while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
178:      KSPLogResidualHistory(ksp,res);
179:      lgmres->it = (loc_it - 1);
180:      KSPMonitor(ksp,ksp->its,res);

182:     /* see if more space is needed for work vectors */
183:     if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
184:        LGMRESGetNewVectors(ksp,loc_it+1);
185:       /* (loc_it+1) is passed in as number of the first vector that should
186:          be allocated */
187:     }

189:     /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
190:     if (loc_it < it_arnoldi) { /* Arnoldi */
191:        KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
192:     } else { /*aug step */
193:        order = loc_it - it_arnoldi + 1; /* which aug step */
194:        for (ii=0; ii<aug_dim; ii++) {
195:            if (lgmres->aug_order[ii] == order) {
196:               spot = ii;
197:               break; /* must have this because there will be duplicates before aug_ct = aug_dim */
198:             }
199:         }

201:        VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
202:        /*note: an alternate implementation choice would be to only save the AUGVECS and
203:          not A_AUGVEC and then apply the PC here to the augvec */
204:     }

206:     /* update hessenberg matrix and do Gram-Schmidt - new direction is in
207:        VEC_VV(1+loc_it)*/
208:     (*lgmres->orthog)(ksp,loc_it);

210:     /* new entry in hessenburg is the 2-norm of our new direction */
211:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
212:     *HH(loc_it+1,loc_it)   = tt;
213:     *HES(loc_it+1,loc_it)  = tt;


216:     /* check for the happy breakdown */
217:     hapbnd  = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration  */
218:     if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
219:     if (tt > hapbnd) {
220:        tmp = 1.0/tt;
221:        VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
222:     } else {
223:        PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
224:        hapend = PETSC_TRUE;
225:     }

227:     /* Now apply rotations to new col of hessenberg (and right side of system), 
228:        calculate new rotation, and get new residual norm at the same time*/
229:     LGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
230:     if (ksp->reason) break;

232:     loc_it++;
233:     lgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
234: 
235:     PetscObjectTakeAccess(ksp);
236:     ksp->its++;
237:     ksp->rnorm = res;
238:     PetscObjectGrantAccess(ksp);

240:     (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);

242:     /* Catch error in happy breakdown and signal convergence and break from loop */
243:     if (hapend) {
244:       if (!ksp->reason) {
245:         SETERRQ1(0,"You reached the happy break down,but convergence was not indicated. Residual norm = %G",res);
246:       }
247:       break;
248:     }
249:   }
250:   /* END OF ITERATION LOOP */

252:   KSPLogResidualHistory(ksp,res);

254:   /* Monitor if we know that we will not return for a restart */
255:   if (ksp->reason || ksp->its >= max_it) {
256:     KSPMonitor(ksp, ksp->its, res);
257:   }

259:   if (itcount) *itcount    = loc_it;

261:   /*
262:     Down here we have to solve for the "best" coefficients of the Krylov
263:     columns, add the solution values together, and possibly unwind the
264:     preconditioning from the solution
265:    */
266: 
267:   /* Form the solution (or the solution so far) */
268:   /* Note: must pass in (loc_it-1) for iteration count so that BuildLgmresSoln
269:      properly navigates */

271:   BuildLgmresSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);


274:   /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
275:      only if we will be restarting (i.e. this cycle performed it_total
276:      iterations)  */
277:   if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {

279:      /*AUG_TEMP contains the new augmentation vector (assigned in  BuildLgmresSoln) */
280:     if (!lgmres->aug_ct) {
281:         spot = 0;
282:         lgmres->aug_ct++;
283:      } else if (lgmres->aug_ct < aug_dim) {
284:         spot = lgmres->aug_ct;
285:         lgmres->aug_ct++;
286:      } else { /* truncate */
287:         for (ii=0; ii<aug_dim; ii++) {
288:            if (lgmres->aug_order[ii] == aug_dim) {
289:               spot = ii;
290:             }
291:         }
292:      }

294: 

296:      VecCopy(AUG_TEMP, AUGVEC(spot));
297:      /*need to normalize */
298:      VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
299:      inv_tmp_norm = 1.0/tmp_norm;
300:      VecScale(AUGVEC(spot),inv_tmp_norm);

302:      /*set new aug vector to order 1  - move all others back one */
303:      for (ii=0; ii < aug_dim; ii++) {
304:         AUG_ORDER(ii)++;
305:      }
306:      AUG_ORDER(spot) = 1;

308:      /*now add the A*aug vector to A_AUGVEC(spot)  - this is independ. of preconditioning type*/
309:      /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */

311: 
312:      /* first do H+*y */
313:      VecSet(AUG_TEMP,0.0);
314:      VecGetArray(AUG_TEMP, &avec);
315:      for (ii=0; ii < it_total + 1; ii++) {
316:         for (jj=0; jj <= ii+1; jj++) {
317:            avec[jj] += *HES(jj ,ii) * *GRS(ii);
318:         }
319:      }

321:      /*now multiply result by V+ */
322:      VecSet(VEC_TEMP,0.0);
323:      VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
324:      VecRestoreArray(AUG_TEMP, &avec);
325: 
326:      /*copy answer to aug location  and scale*/
327:      VecCopy(VEC_TEMP,  A_AUGVEC(spot));
328:      VecScale(A_AUGVEC(spot),inv_tmp_norm);
329:   }
330:   return(0);
331: }

333: /*  
334:     KSPSolve_LGMRES - This routine applies the LGMRES method.


337:    Input Parameter:
338: .     ksp - the Krylov space object that was set to use lgmres

340:    Output Parameter:
341: .     outits - number of iterations used

343: */

347: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
348: {
350:   PetscInt       cycle_its; /* iterations done in a call to LGMREScycle */
351:   PetscInt       itcount;   /* running total of iterations, incl. those in restarts */
352:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
353:   PetscTruth     guess_zero = ksp->guess_zero;
354:   PetscInt       ii;        /*LGMRES_MOD variable */

357:   if (ksp->calc_sings && !lgmres->Rsvd) {
358:      SETERRQ(PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
359:   }
360:   PetscObjectTakeAccess(ksp);
361:   ksp->its        = 0;
362:   lgmres->aug_ct  = 0;
363:   lgmres->matvecs = 0;
364:   PetscObjectGrantAccess(ksp);

366:   /* initialize */
367:   itcount  = 0;
368:   ksp->reason = KSP_CONVERGED_ITERATING;
369:   /*LGMRES_MOD*/
370:   for (ii=0; ii<lgmres->aug_dim; ii++) {
371:      lgmres->aug_order[ii] = 0;
372:   }

374:   while (!ksp->reason) {
375:      /* calc residual - puts in VEC_VV(0) */
376:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
377:     LGMREScycle(&cycle_its,ksp);
378:     itcount += cycle_its;
379:     if (itcount >= ksp->max_it) {
380:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
381:       break;
382:     }
383:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
384:   }
385:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
386:   return(0);
387: }

390: /*

392:    KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.

394: */
397: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
398: {
399:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;

403:   PetscFree(lgmres->augvecs);
404:   if (lgmres->augwork_alloc) {
405:     VecDestroyVecs(lgmres->augvecs_user_work[0],lgmres->augwork_alloc);
406:   }
407:   PetscFree(lgmres->augvecs_user_work);
408:   PetscFree(lgmres->aug_order);
409:   KSPDestroy_GMRES(ksp);
410:   return(0);
411: }

413: /*
414:     BuildLgmresSoln - create the solution from the starting vector and the
415:                       current iterates.

417:     Input parameters:
418:         nrs - work area of size it + 1.
419:         vguess  - index of initial guess
420:         vdest - index of result.  Note that vguess may == vdest (replace
421:                 guess with the solution).
422:         it - HH upper triangular part is a block of size (it+1) x (it+1)  

424:      This is an internal routine that knows about the LGMRES internals.
425:  */
428: static PetscErrorCode BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
429: {
430:   PetscScalar    tt;
432:   PetscInt       ii,k,j;
433:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)(ksp->data);
434:   /*LGMRES_MOD */
435:   PetscInt       it_arnoldi, it_aug;
436:   PetscInt       jj, spot = 0;

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

441:   /* If it is < 0, no lgmres steps have been performed */
442:   if (it < 0) {
443:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
444:     return(0);
445:   }

447:   /* so (it+1) lgmres steps HAVE been performed */

449:   /* LGMRES_MOD - determine if we need to use augvecs for the soln  - do not assume that
450:      this is called after the total its allowed for an approx space */
451:    if (lgmres->approx_constant) {
452:      it_arnoldi = lgmres->max_k - lgmres->aug_ct;
453:    } else {
454:      it_arnoldi = lgmres->max_k - lgmres->aug_dim;
455:    }
456:    if (it_arnoldi >= it +1) {
457:       it_aug = 0;
458:       it_arnoldi = it+1;
459:    } else {
460:       it_aug = (it + 1) - it_arnoldi;
461:    }

463:   /* now it_arnoldi indicates the number of matvecs that took place */
464:   lgmres->matvecs += it_arnoldi;

466: 
467:   /* solve the upper triangular system - GRS is the right side and HH is 
468:      the upper triangular matrix  - put soln in nrs */
469:   if (*HH(it,it) == 0.0) SETERRQ2(PETSC_ERR_CONV_FAILED,"HH(it,it) is identically zero; it = %D GRS(it) = %G",it,PetscAbsScalar(*GRS(it)));
470:   if (*HH(it,it) != 0.0) {
471:      nrs[it] = *GRS(it) / *HH(it,it);
472:   } else {
473:      nrs[it] = 0.0;
474:   }

476:   for (ii=1; ii<=it; ii++) {
477:     k   = it - ii;
478:     tt  = *GRS(k);
479:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
480:     nrs[k]   = tt / *HH(k,k);
481:   }

483:   /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
484:   VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */

486:   /*LGMRES_MOD - if augmenting has happened we need to form the solution 
487:     using the augvecs */
488:   if (!it_aug) { /* all its are from arnoldi */
489:      VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
490:   } else { /*use aug vecs */
491:      /*first do regular krylov directions */
492:      VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
493:      /*now add augmented portions - add contribution of aug vectors one at a time*/


496:      for (ii=0; ii<it_aug; ii++) {
497:         for (jj=0; jj<lgmres->aug_dim; jj++) {
498:            if (lgmres->aug_order[jj] == (ii+1)) {
499:               spot = jj;
500:               break; /* must have this because there will be duplicates before aug_ct = aug_dim */
501:             }
502:         }
503:         VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
504:       }
505:   }
506:   /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
507:      preconditioner is "unwound" from right-precondtioning*/
508:   VecCopy(VEC_TEMP, AUG_TEMP);

510:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);

512:   /* add solution to previous solution */
513:   /* put updated solution into vdest.*/
514:   VecCopy(vguess,vdest);
515:   VecAXPY(vdest,1.0,VEC_TEMP);

517:   return(0);
518: }

520: /*

522:     LGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
523:                             Return new residual.

525:     input parameters:

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

532:     output parameters:
533: .        res - the new residual
534:         
535:  */
538: static PetscErrorCode LGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
539: {
540:   PetscScalar   *hh,*cc,*ss,tt;
541:   PetscInt      j;
542:   KSP_LGMRES    *lgmres = (KSP_LGMRES *)(ksp->data);

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

550:   /* Apply all the previously computed plane rotations to the new column
551:      of the Hessenberg matrix */
552:   /* Note: this uses the rotation [conj(c)  s ; -s   c], c= cos(theta), s= sin(theta) */

554:   for (j=1; j<=it; j++) {
555:     tt  = *hh;
556: #if defined(PETSC_USE_COMPLEX)
557:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
558: #else
559:     *hh = *cc * tt + *ss * *(hh+1);
560: #endif
561:     hh++;
562:     *hh = *cc++ * *hh - (*ss++ * tt);
563:     /* hh, cc, and ss have all been incremented one by end of loop */
564:   }

566:   /*
567:     compute the new plane rotation, and apply it to:
568:      1) the right-hand-side of the Hessenberg system (GRS)
569:         note: it affects GRS(it) and GRS(it+1)
570:      2) the new column of the Hessenberg matrix
571:         note: it affects HH(it,it) which is currently pointed to 
572:         by hh and HH(it+1, it) (*(hh+1))  
573:     thus obtaining the updated value of the residual...
574:   */

576:   /* compute new plane rotation */

578:   if (!hapend) {
579: #if defined(PETSC_USE_COMPLEX)
580:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
581: #else
582:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
583: #endif
584:     if (tt == 0.0) {
585:       ksp->reason = KSP_DIVERGED_NULL;
586:       return(0);
587:     }
588:     *cc       = *hh / tt;   /* new cosine value */
589:     *ss       = *(hh+1) / tt;  /* new sine value */

591:     /* apply to 1) and 2) */
592:     *GRS(it+1) = - (*ss * *GRS(it));
593: #if defined(PETSC_USE_COMPLEX)
594:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
595:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
596: #else
597:     *GRS(it)   = *cc * *GRS(it);
598:     *hh        = *cc * *hh + *ss * *(hh+1);
599: #endif

601:     /* residual is the last element (it+1) of right-hand side! */
602:     *res      = PetscAbsScalar(*GRS(it+1));

604:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
605:             another rotation matrix (so RH doesn't change).  The new residual is 
606:             always the new sine term times the residual from last time (GRS(it)), 
607:             but now the new sine rotation would be zero...so the residual should
608:             be zero...so we will multiply "zero" by the last residual.  This might
609:             not be exactly what we want to do here -could just return "zero". */
610: 
611:     *res = 0.0;
612:   }
613:   return(0);
614: }

616: /*

618:    LGMRESGetNewVectors - This routine allocates more work vectors, starting from 
619:                          VEC_VV(it) 
620:                          
621: */
624: static PetscErrorCode LGMRESGetNewVectors(KSP ksp,PetscInt it)
625: {
626:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
627:   PetscInt       nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
628:   PetscInt       nalloc;                      /* number to allocate */
630:   PetscInt       k;
631: 
633:   nalloc = lgmres->delta_allocate; /* number of vectors to allocate 
634:                                       in a single chunk */

636:   /* Adjust the number to allocate to make sure that we don't exceed the
637:      number of available slots (lgmres->vecs_allocated)*/
638:   if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
639:     nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
640:   }
641:   if (!nalloc) return(0);

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

645:   /* work vectors */
646:   KSPGetVecs(ksp,nalloc,&lgmres->user_work[nwork],0,PETSC_NULL);
647:   PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
648:   /* specify size of chunk allocated */
649:   lgmres->mwork_alloc[nwork] = nalloc;

651:   for (k=0; k < nalloc; k++) {
652:     lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
653:   }
654: 

656:   /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
657: 

659:   /* increment the number of work vector chunks */
660:   lgmres->nwork_alloc++;
661:   return(0);
662: }

664: /* 

666:    KSPBuildSolution_LGMRES

668:      Input Parameter:
669: .     ksp - the Krylov space object
670: .     ptr-

672:    Output Parameter:
673: .     result - the solution

675:    Note: this calls BuildLgmresSoln - the same function that LGMREScycle
676:    calls directly.  

678: */
681: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
682: {
683:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;

687:   if (!ptr) {
688:     if (!lgmres->sol_temp) {
689:       VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
690:       PetscLogObjectParent(ksp,lgmres->sol_temp);
691:     }
692:     ptr = lgmres->sol_temp;
693:   }
694:   if (!lgmres->nrs) {
695:     /* allocate the work area */
696:     PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
697:     PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
698:   }
699: 
700:   BuildLgmresSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
701:   if (result) *result = ptr;
702: 
703:   return(0);
704: }


710: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
711: {
712:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
714:   PetscTruth     iascii;

717:   KSPView_GMRES(ksp,viewer);
718:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
719:   if (iascii) {
720:     /*LGMRES_MOD */
721:     PetscViewerASCIIPrintf(viewer,"  LGMRES: aug. dimension=%D\n",lgmres->aug_dim);
722:     if (lgmres->approx_constant) {
723:        PetscViewerASCIIPrintf(viewer,"  LGMRES: approx. space size was kept constant.\n");
724:     }
725:     PetscViewerASCIIPrintf(viewer,"  LGMRES: number of matvecs=%D\n",lgmres->matvecs);
726:   } else {
727:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
728:   }
729:   return(0);
730: }


736: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp)
737: {
739:   PetscInt       aug;
740:   KSP_LGMRES     *lgmres = (KSP_LGMRES*) ksp->data;
741:   PetscTruth     flg = PETSC_FALSE;

744:   KSPSetFromOptions_GMRES(ksp);
745:   PetscOptionsHead("KSP LGMRES Options");
746:   PetscOptionsTruth("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",flg,&flg,PETSC_NULL);
747:     if (flg) { lgmres->approx_constant = 1; }
748:     PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
749:     if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
750:   PetscOptionsTail();
751:   return(0);
752: }


755: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
756: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);

758: /*functions for extra lgmres options here*/
762: PetscErrorCode  KSPLGMRESSetConstant_LGMRES(KSP ksp)
763: {
764:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
766:   lgmres->approx_constant = 1;
767:   return(0);
768: }

774: PetscErrorCode  KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
775: {
776:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

779:   if (aug_dim < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
780:   if (aug_dim > (lgmres->max_k -1))  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
781:   lgmres->aug_dim = aug_dim;
782:   return(0);
783: }


787: /* end new lgmres functions */


790: /* use these options from gmres */
792: EXTERN PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP,double);
793: EXTERN PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP);
794: EXTERN PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP,PetscInt);
795: EXTERN PetscErrorCode  KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));
796: EXTERN PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);

799: /*MC
800:     KSPLGMRES - Augments the standard GMRES approximation space with approximations to
801:                 the error from previous restart cycles.

803:   Options Database Keys:
804: +   -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
805: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
806: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
807:                             vectors are allocated as needed)
808: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
809: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
810: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
811:                                   stability of the classical Gram-Schmidt  orthogonalization.
812: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
813: .   -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
814: -   -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)

816:    Described in:
817:     A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
818:     accelerating the convergence of restarted GMRES. SIAM
819:     Journal on Matrix Analysis and Applications, 26 (2005), pp. 962-984.

821:     To run LGMRES(m, k) as described in the above paper, use:
822:        -ksp_gmres_restart <m+k>
823:        -ksp_lgmres_augment <k>

825:   Level: beginner

827:   Notes:  This object is subclassed off of KSPGMRES

829:   Contributed by: Allison Baker

831: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
832:           KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
833:           KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
834:           KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
835:           KSPGMRESSetConstant()

837: M*/

842: PetscErrorCode  KSPCreate_LGMRES(KSP ksp)
843: {
844:   KSP_LGMRES     *lgmres;

848:   PetscNewLog(ksp,KSP_LGMRES,&lgmres);
849:   ksp->data                              = (void*)lgmres;
850:   ksp->ops->buildsolution                = KSPBuildSolution_LGMRES;

852:   ksp->ops->setup                        = KSPSetUp_LGMRES;
853:   ksp->ops->solve                        = KSPSolve_LGMRES;
854:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
855:   ksp->ops->view                         = KSPView_LGMRES;
856:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
857:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
858:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

860:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
861:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
862:                                      KSPGMRESSetPreAllocateVectors_GMRES);
863:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
864:                                     "KSPGMRESSetOrthogonalization_GMRES",
865:                                      KSPGMRESSetOrthogonalization_GMRES);
866:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
867:                                     "KSPGMRESSetRestart_GMRES",
868:                                      KSPGMRESSetRestart_GMRES);
869:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
870:                                     "KSPGMRESSetHapTol_GMRES",
871:                                      KSPGMRESSetHapTol_GMRES);
872:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
873:                                     "KSPGMRESSetCGSRefinementType_GMRES",
874:                                      KSPGMRESSetCGSRefinementType_GMRES);

876:   /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
877:    PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
878:                                      "KSPLGMRESSetConstant_LGMRES",
879:                                       KSPLGMRESSetConstant_LGMRES);

881:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
882:                                     "KSPLGMRESSetAugDim_LGMRES",
883:                                      KSPLGMRESSetAugDim_LGMRES);
884: 

886:   /*defaults */
887:   lgmres->haptol              = 1.0e-30;
888:   lgmres->q_preallocate       = 0;
889:   lgmres->delta_allocate      = LGMRES_DELTA_DIRECTIONS;
890:   lgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
891:   lgmres->nrs                 = 0;
892:   lgmres->sol_temp            = 0;
893:   lgmres->max_k               = LGMRES_DEFAULT_MAXK;
894:   lgmres->Rsvd                = 0;
895:   lgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
896:   lgmres->orthogwork          = 0;
897:   /*LGMRES_MOD - new defaults */
898:   lgmres->aug_dim             = LGMRES_DEFAULT_AUGDIM;
899:   lgmres->aug_ct              = 0; /* start with no aug vectors */
900:   lgmres->approx_constant     = 0;
901:   lgmres->matvecs             = 0;

903:   return(0);
904: }