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:   PetscMalloc((max_k+1)*sizeof(PetscScalar),&lgmres->hwork);
 71:   PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
 72:   for (k=0; k<lgmres->aug_vv_allocated; k++) {
 73:     lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
 74:   }
 75:   return(0);
 76: }


 79: /*

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

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

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

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


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

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

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

140:   it_total =  it_arnoldi + lgmres->aug_ct;

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

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

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

161:   ksp->rnorm = res;


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

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


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

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

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

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

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

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


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

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

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

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

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

253:   KSPLogResidualHistory(ksp,res);

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

260:   if (itcount) *itcount    = loc_it;

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

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


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

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

295: 

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

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

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

312: 
313:      /* first do H+*y */
314:      avec = lgmres->hwork;
315:      PetscMemzero(avec,(it_total+1)*sizeof(*avec));
316:      for (ii=0; ii < it_total + 1; ii++) {
317:         for (jj=0; jj <= ii+1 && jj < it_total+1; jj++) {
318:            avec[jj] += *HES(jj ,ii) * *GRS(ii);
319:         }
320:      }

322:      /*now multiply result by V+ */
323:      VecSet(VEC_TEMP,0.0);
324:      VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/

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:   PetscFree(lgmres->hwork);
410:   KSPDestroy_GMRES(ksp);
411:   return(0);
412: }

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

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

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

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

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

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

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

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

467: 
468:   /* solve the upper triangular system - GRS is the right side and HH is 
469:      the upper triangular matrix  - put soln in nrs */
470:   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)));
471:   if (*HH(it,it) != 0.0) {
472:      nrs[it] = *GRS(it) / *HH(it,it);
473:   } else {
474:      nrs[it] = 0.0;
475:   }

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

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

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


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

511:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);

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

518:   return(0);
519: }

521: /*

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

526:     input parameters:

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

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

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

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

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

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

577:   /* compute new plane rotation */

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

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

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

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

617: /*

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

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

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

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

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

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

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

665: /* 

667:    KSPBuildSolution_LGMRES

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

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

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

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

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


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

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


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

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


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

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

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

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


788: /* end new lgmres functions */


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

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

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

817:     To run LGMRES(m, k) as described in the above paper, use:
818:        -ksp_gmres_restart <m+k>
819:        -ksp_lgmres_augment <k>

821:   Level: beginner

823:    Notes: Supports both left and right preconditioning, but not symmetric.

825:    References:
826:     A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
827:     accelerating the convergence of restarted GMRES. SIAM
828:     Journal on Matrix Analysis and Applications, 26 (2005), pp. 962-984.

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

832:   Contributed by: Allison Baker

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

840: M*/

845: PetscErrorCode  KSPCreate_LGMRES(KSP ksp)
846: {
847:   KSP_LGMRES     *lgmres;

851:   PetscNewLog(ksp,KSP_LGMRES,&lgmres);
852:   ksp->data                              = (void*)lgmres;
853:   ksp->ops->buildsolution                = KSPBuildSolution_LGMRES;

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

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

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

884:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
885:                                     "KSPLGMRESSetAugDim_LGMRES",
886:                                      KSPLGMRESSetAugDim_LGMRES);
887: 

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

906:   return(0);
907: }