Actual source code: lgmres.c

  1: #define PETSCKSP_DLL

 3:  #include ../src/ksp/ksp/impls/gmres/lgmres/lgmresp.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: }

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

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

 40: */
 43: PetscErrorCode    KSPSetUp_LGMRES(KSP ksp)
 44: {
 45:   PetscInt       size,hh,hes,rs,cc;
 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:   hh            = (max_k + 2) * (max_k + 1);
 57:   hes           = (max_k + 1) * (max_k + 1);
 58:   rs            = (max_k + 2);
 59:   cc            = (max_k + 1);  /* SS and CC are the same size */
 60:   size          = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);

 62:   /* Allocate space and set pointers to beginning */
 63:   PetscMalloc(size,&lgmres->hh_origin);
 64:   PetscMemzero(lgmres->hh_origin,size);
 65:   PetscLogObjectMemory(ksp,size);  /* HH - modified (by plane rotations) hessenburg */
 66:   lgmres->hes_origin = lgmres->hh_origin + hh;     /* HES - unmodified hessenburg */
 67:   lgmres->rs_origin  = lgmres->hes_origin + hes;   /* RS - the right-hand-side of the 
 68:                                                       Hessenberg system */
 69:   lgmres->cc_origin  = lgmres->rs_origin + rs;     /* CC - cosines for rotations */
 70:   lgmres->ss_origin  = lgmres->cc_origin + cc;     /* SS - sines for rotations */

 72:   if (ksp->calc_sings) {
 73:     /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
 74:     size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
 75:     PetscMalloc(size,&lgmres->Rsvd);
 76:     PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&lgmres->Dsvd);
 77:     PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
 78:   }

 80:   /* Allocate array to hold pointers to user vectors.  Note that we need
 81:   we need it+1 vectors, and it <= max_k)  - vec_offset indicates some initial work vectors*/
 82:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&lgmres->vecs);
 83:   lgmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
 84:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&lgmres->user_work);
 85:   PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&lgmres->mwork_alloc);
 86:   PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));

 88:   /* LGMRES_MOD: need array of pointers to augvecs*/
 89:   PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs);
 90:   lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
 91:   PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs_user_work);
 92:   PetscMalloc(aug_dim*sizeof(PetscInt),&lgmres->aug_order);
 93:   PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));

 95: 
 96:  /* if q_preallocate = 0 then only allocate one "chunk" of space (for 
 97:      5 vectors) - additional will then be allocated from LGMREScycle() 
 98:      as needed.  Otherwise, allocate all of the space that could be needed */
 99:   if (lgmres->q_preallocate) {
100:     lgmres->vv_allocated   = VEC_OFFSET + 2 + max_k;
101:     KSPGetVecs(ksp,lgmres->vv_allocated,&lgmres->user_work[0],0,PETSC_NULL);
102:     PetscLogObjectParents(ksp,lgmres->vv_allocated,lgmres->user_work[0]);
103:     lgmres->mwork_alloc[0] = lgmres->vv_allocated;
104:     lgmres->nwork_alloc    = 1;
105:     for (k=0; k<lgmres->vv_allocated; k++) {
106:       lgmres->vecs[k] = lgmres->user_work[0][k];
107:     }
108:   } else {
109:     lgmres->vv_allocated    = 5;
110:     KSPGetVecs(ksp,5,&lgmres->user_work[0],0,PETSC_NULL);
111:     PetscLogObjectParents(ksp,5,lgmres->user_work[0]);
112:     lgmres->mwork_alloc[0]  = 5;
113:     lgmres->nwork_alloc     = 1;
114:     for (k=0; k<lgmres->vv_allocated; k++) {
115:       lgmres->vecs[k] = lgmres->user_work[0][k];
116:     }
117:   }
118:   /* LGMRES_MOD - for now we will preallocate the augvecs - because aug_dim << restart
119:      ... also keep in mind that we need to keep augvecs from cycle to cycle*/
120:   lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
121:   lgmres->augwork_alloc =  2* aug_dim + AUG_OFFSET;
122:   KSPGetVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,PETSC_NULL);
123:   PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
124:   for (k=0; k<lgmres->aug_vv_allocated; k++) {
125:     lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
126:   }
127:   return(0);
128: }


131: /*

133:     LGMRESCycle - Run lgmres, possibly with restart.  Return residual 
134:                   history if requested.

136:     input parameters:
137: .         lgmres  - structure containing parameters and work areas

139:     output parameters:
140: .        nres    - residuals (from preconditioned system) at each step.
141:                   If restarting, consider passing nres+it.  If null, 
142:                   ignored
143: .        itcount - number of iterations used.   nres[0] to nres[itcount]
144:                   are defined.  If null, ignored.  If null, ignored.
145: .        converged - 0 if not converged

147:                   
148:     Notes:
149:     On entry, the value in vector VEC_VV(0) should be 
150:     the initial residual.


153:  */
156: PetscErrorCode LGMREScycle(PetscInt *itcount,KSP ksp)
157: {

159:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)(ksp->data);
160:   PetscReal      res_norm, res;
161:   PetscReal      hapbnd, tt;
162:   PetscScalar    tmp;
163:   PetscTruth     hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
165:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
166:   PetscInt       max_k = lgmres->max_k; /* max approx space size */
167:   PetscInt       max_it = ksp->max_it;  /* max # of overall iterations for the method */
168:   /* LGMRES_MOD - new variables*/
169:   PetscInt       aug_dim = lgmres->aug_dim;
170:   PetscInt       spot = 0;
171:   PetscInt       order = 0;
172:   PetscInt       it_arnoldi;             /* number of arnoldi steps to take */
173:   PetscInt       it_total;               /* total number of its to take (=approx space size)*/
174:   PetscInt       ii, jj;
175:   PetscReal      tmp_norm;
176:   PetscScalar    inv_tmp_norm;
177:   PetscScalar    *avec;

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

184:   /* LGMRES_MOD: determine number of arnoldi steps to take */
185:   /* if approx_constant then we keep the space the same size even if 
186:      we don't have the full number of aug vectors yet*/
187:   if (lgmres->approx_constant) {
188:      it_arnoldi = max_k - lgmres->aug_ct;
189:   } else {
190:       it_arnoldi = max_k - aug_dim;
191:   }

193:   it_total =  it_arnoldi + lgmres->aug_ct;

195:   /* initial residual is in VEC_VV(0)  - compute its norm*/
196:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
197:   res    = res_norm;
198: 
199:   /* first entry in right-hand-side of hessenberg system is just 
200:      the initial residual norm */
201:   *GRS(0) = res_norm;

203:  /* check for the convergence */
204:   if (!res) {
205:      if (itcount) *itcount = 0;
206:      ksp->reason = KSP_CONVERGED_ATOL;
207:      PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
208:      return(0);
209:   }

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

214:   ksp->rnorm = res;


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

223: 
224:   /* MAIN ITERATION LOOP BEGINNING*/


227:   /* keep iterating until we have converged OR generated the max number
228:      of directions OR reached the max number of iterations for the method */
229:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
230: 
231:   while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
232:      KSPLogResidualHistory(ksp,res);
233:      lgmres->it = (loc_it - 1);
234:      KSPMonitor(ksp,ksp->its,res);

236:     /* see if more space is needed for work vectors */
237:     if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
238:        LGMRESGetNewVectors(ksp,loc_it+1);
239:       /* (loc_it+1) is passed in as number of the first vector that should
240:          be allocated */
241:     }

243:     /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
244:     if (loc_it < it_arnoldi) { /* Arnoldi */
245:        KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
246:     } else { /*aug step */
247:        order = loc_it - it_arnoldi + 1; /* which aug step */
248:        for (ii=0; ii<aug_dim; ii++) {
249:            if (lgmres->aug_order[ii] == order) {
250:               spot = ii;
251:               break; /* must have this because there will be duplicates before aug_ct = aug_dim */
252:             }
253:         }

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

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

264:     /* new entry in hessenburg is the 2-norm of our new direction */
265:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
266:     *HH(loc_it+1,loc_it)   = tt;
267:     *HES(loc_it+1,loc_it)  = tt;


270:     /* check for the happy breakdown */
271:     hapbnd  = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration  */
272:     if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
273:     if (tt > hapbnd) {
274:        tmp = 1.0/tt;
275:        VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
276:     } else {
277:        PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
278:        hapend = PETSC_TRUE;
279:     }

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

286:     loc_it++;
287:     lgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
288: 
289:     PetscObjectTakeAccess(ksp);
290:     ksp->its++;
291:     ksp->rnorm = res;
292:     PetscObjectGrantAccess(ksp);

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

296:     /* Catch error in happy breakdown and signal convergence and break from loop */
297:     if (hapend) {
298:       if (!ksp->reason) {
299:         SETERRQ1(0,"You reached the happy break down,but convergence was not indicated. Residual norm = %G",res);
300:       }
301:       break;
302:     }
303:   }
304:   /* END OF ITERATION LOOP */

306:   KSPLogResidualHistory(ksp,res);

308:   /* Monitor if we know that we will not return for a restart */
309:   if (ksp->reason || ksp->its >= max_it) {
310:     KSPMonitor(ksp, ksp->its, res);
311:   }

313:   if (itcount) *itcount    = loc_it;

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

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


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

333:      /*AUG_TEMP contains the new augmentation vector (assigned in  BuildLgmresSoln) */
334:     if (!lgmres->aug_ct) {
335:         spot = 0;
336:         lgmres->aug_ct++;
337:      } else if (lgmres->aug_ct < aug_dim) {
338:         spot = lgmres->aug_ct;
339:         lgmres->aug_ct++;
340:      } else { /* truncate */
341:         for (ii=0; ii<aug_dim; ii++) {
342:            if (lgmres->aug_order[ii] == aug_dim) {
343:               spot = ii;
344:             }
345:         }
346:      }

348: 

350:      VecCopy(AUG_TEMP, AUGVEC(spot));
351:      /*need to normalize */
352:      VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
353:      inv_tmp_norm = 1.0/tmp_norm;
354:      VecScale(AUGVEC(spot),inv_tmp_norm);

356:      /*set new aug vector to order 1  - move all others back one */
357:      for (ii=0; ii < aug_dim; ii++) {
358:         AUG_ORDER(ii)++;
359:      }
360:      AUG_ORDER(spot) = 1;

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

365: 
366:      /* first do H+*y */
367:      VecSet(AUG_TEMP,0.0);
368:      VecGetArray(AUG_TEMP, &avec);
369:      for (ii=0; ii < it_total + 1; ii++) {
370:         for (jj=0; jj <= ii+1; jj++) {
371:            avec[jj] += *HES(jj ,ii) * *GRS(ii);
372:         }
373:      }

375:      /*now multiply result by V+ */
376:      VecSet(VEC_TEMP,0.0);
377:      VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
378:      VecRestoreArray(AUG_TEMP, &avec);
379: 
380:      /*copy answer to aug location  and scale*/
381:      VecCopy(VEC_TEMP,  A_AUGVEC(spot));
382:      VecScale(A_AUGVEC(spot),inv_tmp_norm);


385:   }
386:   return(0);
387: }

389: /*  
390:     KSPSolve_LGMRES - This routine applies the LGMRES method.


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

396:    Output Parameter:
397: .     outits - number of iterations used

399: */

403: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
404: {
406:   PetscInt       cycle_its; /* iterations done in a call to LGMREScycle */
407:   PetscInt       itcount;   /* running total of iterations, incl. those in restarts */
408:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
409:   PetscTruth     guess_zero = ksp->guess_zero;
410:   PetscInt       ii;        /*LGMRES_MOD variable */

413:   if (ksp->calc_sings && !lgmres->Rsvd) {
414:      SETERRQ(PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
415:   }
416:   PetscObjectTakeAccess(ksp);
417:   ksp->its        = 0;
418:   lgmres->aug_ct  = 0;
419:   lgmres->matvecs = 0;
420:   PetscObjectGrantAccess(ksp);

422:   /* initialize */
423:   itcount  = 0;
424:   ksp->reason = KSP_CONVERGED_ITERATING;
425:   /*LGMRES_MOD*/
426:   for (ii=0; ii<lgmres->aug_dim; ii++) {
427:      lgmres->aug_order[ii] = 0;
428:   }

430:   while (!ksp->reason) {
431:      /* calc residual - puts in VEC_VV(0) */
432:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
433:     LGMREScycle(&cycle_its,ksp);
434:     itcount += cycle_its;
435:     if (itcount >= ksp->max_it) {
436:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
437:       break;
438:     }
439:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
440:   }
441:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
442:   return(0);
443: }

445: /*

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

449: */
452: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
453: {
454:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;
456:   PetscInt       i;

459:   /* Free the Hessenberg matrices */
460:   PetscFree(lgmres->hh_origin);

462:   /* Free pointers to user variables */
463:   PetscFree(lgmres->vecs);

465:   /*LGMRES_MOD - free pointers for extra vectors */
466:   PetscFree(lgmres->augvecs);

468:   /* free work vectors */
469:   for (i=0; i < lgmres->nwork_alloc; i++) {
470:     VecDestroyVecs(lgmres->user_work[i],lgmres->mwork_alloc[i]);
471:   }
472:   PetscFree(lgmres->user_work);

474:   /*LGMRES_MOD - free aug work vectors also */
475:   /*this was all allocated as one "chunk" */
476:   if (lgmres->augwork_alloc) {
477:     VecDestroyVecs(lgmres->augvecs_user_work[0],lgmres->augwork_alloc);
478:   }
479:   PetscFree(lgmres->augvecs_user_work);
480:   PetscFree(lgmres->aug_order);
481:   PetscFree(lgmres->mwork_alloc);
482:   PetscFree(lgmres->nrs);
483:   if (lgmres->sol_temp) {VecDestroy(lgmres->sol_temp);}
484:   PetscFree(lgmres->Rsvd);
485:   PetscFree(lgmres->Dsvd);
486:   PetscFree(lgmres->orthogwork);
487:   PetscFree(ksp->data);
488:   /* clear composed functions */
489:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C","",PETSC_NULL);
490:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C","",PETSC_NULL);
491:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C","",PETSC_NULL);
492:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C","",PETSC_NULL);
493:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C","",PETSC_NULL);
494:   return(0);
495: }

497: /*
498:     BuildLgmresSoln - create the solution from the starting vector and the
499:                       current iterates.

501:     Input parameters:
502:         nrs - work area of size it + 1.
503:         vguess  - index of initial guess
504:         vdest - index of result.  Note that vguess may == vdest (replace
505:                 guess with the solution).
506:         it - HH upper triangular part is a block of size (it+1) x (it+1)  

508:      This is an internal routine that knows about the LGMRES internals.
509:  */
512: static PetscErrorCode BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
513: {
514:   PetscScalar    tt;
516:   PetscInt       ii,k,j;
517:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)(ksp->data);
518:   /*LGMRES_MOD */
519:   PetscInt       it_arnoldi, it_aug;
520:   PetscInt       jj, spot = 0;

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

525:   /* If it is < 0, no lgmres steps have been performed */
526:   if (it < 0) {
527:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
528:     return(0);
529:   }

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

533:   /* LGMRES_MOD - determine if we need to use augvecs for the soln  - do not assume that
534:      this is called after the total its allowed for an approx space */
535:    if (lgmres->approx_constant) {
536:      it_arnoldi = lgmres->max_k - lgmres->aug_ct;
537:    } else {
538:      it_arnoldi = lgmres->max_k - lgmres->aug_dim;
539:    }
540:    if (it_arnoldi >= it +1) {
541:       it_aug = 0;
542:       it_arnoldi = it+1;
543:    } else {
544:       it_aug = (it + 1) - it_arnoldi;
545:    }

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

550: 
551:   /* solve the upper triangular system - GRS is the right side and HH is 
552:      the upper triangular matrix  - put soln in nrs */
553:   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)));
554:   if (*HH(it,it) != 0.0) {
555:      nrs[it] = *GRS(it) / *HH(it,it);
556:   } else {
557:      nrs[it] = 0.0;
558:   }

560:   for (ii=1; ii<=it; ii++) {
561:     k   = it - ii;
562:     tt  = *GRS(k);
563:     for (j=k+1; j<=it; j++) tt  = tt - *HH(k,j) * nrs[j];
564:     nrs[k]   = tt / *HH(k,k);
565:   }

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

570:   /*LGMRES_MOD - if augmenting has happened we need to form the solution 
571:     using the augvecs */
572:   if (!it_aug) { /* all its are from arnoldi */
573:      VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
574:   } else { /*use aug vecs */
575:      /*first do regular krylov directions */
576:      VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
577:      /*now add augmented portions - add contribution of aug vectors one at a time*/


580:      for (ii=0; ii<it_aug; ii++) {
581:         for (jj=0; jj<lgmres->aug_dim; jj++) {
582:            if (lgmres->aug_order[jj] == (ii+1)) {
583:               spot = jj;
584:               break; /* must have this because there will be duplicates before aug_ct = aug_dim */
585:             }
586:         }
587:         VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
588:       }
589:   }
590:   /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
591:      preconditioner is "unwound" from right-precondtioning*/
592:   VecCopy(VEC_TEMP, AUG_TEMP);

594:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);

596:   /* add solution to previous solution */
597:   /* put updated solution into vdest.*/
598:   if (vdest != vguess) {
599:     VecCopy(VEC_TEMP,vdest);
600:   }
601:   VecAXPY(vdest,1.0,VEC_TEMP);

603:   return(0);
604: }

606: /*

608:     LGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
609:                             Return new residual.

611:     input parameters:

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

618:     output parameters:
619: .        res - the new residual
620:         
621:  */
624: static PetscErrorCode LGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
625: {
626:   PetscScalar   *hh,*cc,*ss,tt;
627:   PetscInt      j;
628:   KSP_LGMRES    *lgmres = (KSP_LGMRES *)(ksp->data);

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

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

640:   for (j=1; j<=it; j++) {
641:     tt  = *hh;
642: #if defined(PETSC_USE_COMPLEX)
643:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
644: #else
645:     *hh = *cc * tt + *ss * *(hh+1);
646: #endif
647:     hh++;
648:     *hh = *cc++ * *hh - (*ss++ * tt);
649:     /* hh, cc, and ss have all been incremented one by end of loop */
650:   }

652:   /*
653:     compute the new plane rotation, and apply it to:
654:      1) the right-hand-side of the Hessenberg system (GRS)
655:         note: it affects GRS(it) and GRS(it+1)
656:      2) the new column of the Hessenberg matrix
657:         note: it affects HH(it,it) which is currently pointed to 
658:         by hh and HH(it+1, it) (*(hh+1))  
659:     thus obtaining the updated value of the residual...
660:   */

662:   /* compute new plane rotation */

664:   if (!hapend) {
665: #if defined(PETSC_USE_COMPLEX)
666:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
667: #else
668:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
669: #endif
670:     if (tt == 0.0) {
671:       ksp->reason = KSP_DIVERGED_NULL;
672:       return(0);
673:     }
674:     *cc       = *hh / tt;   /* new cosine value */
675:     *ss       = *(hh+1) / tt;  /* new sine value */

677:     /* apply to 1) and 2) */
678:     *GRS(it+1) = - (*ss * *GRS(it));
679: #if defined(PETSC_USE_COMPLEX)
680:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
681:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
682: #else
683:     *GRS(it)   = *cc * *GRS(it);
684:     *hh        = *cc * *hh + *ss * *(hh+1);
685: #endif

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

690:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply 
691:             another rotation matrix (so RH doesn't change).  The new residual is 
692:             always the new sine term times the residual from last time (GRS(it)), 
693:             but now the new sine rotation would be zero...so the residual should
694:             be zero...so we will multiply "zero" by the last residual.  This might
695:             not be exactly what we want to do here -could just return "zero". */
696: 
697:     *res = 0.0;
698:   }
699:   return(0);
700: }

702: /*

704:    LGMRESGetNewVectors - This routine allocates more work vectors, starting from 
705:                          VEC_VV(it) 
706:                          
707: */
710: static PetscErrorCode LGMRESGetNewVectors(KSP ksp,PetscInt it)
711: {
712:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
713:   PetscInt       nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
714:   PetscInt       nalloc;                      /* number to allocate */
716:   PetscInt       k;
717: 
719:   nalloc = lgmres->delta_allocate; /* number of vectors to allocate 
720:                                       in a single chunk */

722:   /* Adjust the number to allocate to make sure that we don't exceed the
723:      number of available slots (lgmres->vecs_allocated)*/
724:   if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
725:     nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
726:   }
727:   if (!nalloc) return(0);

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

731:   /* work vectors */
732:   KSPGetVecs(ksp,nalloc,&lgmres->user_work[nwork],0,PETSC_NULL);
733:   PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
734:   /* specify size of chunk allocated */
735:   lgmres->mwork_alloc[nwork] = nalloc;

737:   for (k=0; k < nalloc; k++) {
738:     lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
739:   }
740: 

742:   /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
743: 

745:   /* increment the number of work vector chunks */
746:   lgmres->nwork_alloc++;
747:   return(0);
748: }

750: /* 

752:    KSPBuildSolution_LGMRES

754:      Input Parameter:
755: .     ksp - the Krylov space object
756: .     ptr-

758:    Output Parameter:
759: .     result - the solution

761:    Note: this calls BuildLgmresSoln - the same function that LGMREScycle
762:    calls directly.  

764: */
767: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
768: {
769:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;

773:   if (!ptr) {
774:     if (!lgmres->sol_temp) {
775:       VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
776:       PetscLogObjectParent(ksp,lgmres->sol_temp);
777:     }
778:     ptr = lgmres->sol_temp;
779:   }
780:   if (!lgmres->nrs) {
781:     /* allocate the work area */
782:     PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
783:     PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
784:   }
785: 
786:   BuildLgmresSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
787:   if (result) *result = ptr;
788: 
789:   return(0);
790: }


796: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
797: {
798:   KSP_LGMRES     *lgmres = (KSP_LGMRES *)ksp->data;
800:   PetscTruth     iascii;

803:   KSPView_GMRES(ksp,viewer);
804:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
805:   if (iascii) {
806:     /*LGMRES_MOD */
807:     PetscViewerASCIIPrintf(viewer,"  LGMRES: aug. dimension=%D\n",lgmres->aug_dim);
808:     if (lgmres->approx_constant) {
809:        PetscViewerASCIIPrintf(viewer,"  LGMRES: approx. space size was kept constant.\n");
810:     }
811:     PetscViewerASCIIPrintf(viewer,"  LGMRES: number of matvecs=%D\n",lgmres->matvecs);
812:   } else {
813:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
814:   }
815:   return(0);
816: }


822: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp)
823: {
825:   PetscInt       aug;
826:   KSP_LGMRES     *lgmres = (KSP_LGMRES*) ksp->data;
827:   PetscTruth     flg;

830:   KSPSetFromOptions_GMRES(ksp);
831:   PetscOptionsHead("KSP LGMRES Options");
832:      PetscOptionsName("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",&flg);
833:     if (flg) { lgmres->approx_constant = 1; }
834:     PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
835:     if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
836:   PetscOptionsTail();
837:   return(0);
838: }


841: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
842: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);

844: /*functions for extra lgmres options here*/
848: PetscErrorCode  KSPLGMRESSetConstant_LGMRES(KSP ksp)
849: {
850:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
852:   lgmres->approx_constant = 1;
853:   return(0);
854: }

860: PetscErrorCode  KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
861: {
862:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

865:   if (aug_dim < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
866:   if (aug_dim > (lgmres->max_k -1))  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
867:   lgmres->aug_dim = aug_dim;
868:   return(0);
869: }


873: /* end new lgmres functions */


876: /* use these options from gmres */
878: EXTERN PetscErrorCode  KSPGMRESSetHapTol_GMRES(KSP,double);
879: EXTERN PetscErrorCode  KSPGMRESSetPreAllocateVectors_GMRES(KSP);
880: EXTERN PetscErrorCode  KSPGMRESSetRestart_GMRES(KSP,PetscInt);
881: EXTERN PetscErrorCode  KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));
882: EXTERN PetscErrorCode  KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);

885: /*MC
886:     KSPLGMRES - Augments the standard GMRES approximation space with approximations to
887:                 the error from previous restart cycles.

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

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

907:     To run LGMRES(m, k) as described in the above paper, use:
908:        -ksp_gmres_restart <m+k>
909:        -ksp_lgmres_augment <k>

911:   Level: beginner

913:   Notes:  This object is subclassed off of KSPGMRES

915:   Contributed by: Allison Baker

917: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
918:           KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
919:           KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
920:           KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
921:           KSPGMRESSetConstant()

923: M*/

928: PetscErrorCode  KSPCreate_LGMRES(KSP ksp)
929: {
930:   KSP_LGMRES     *lgmres;

934:   PetscNewLog(ksp,KSP_LGMRES,&lgmres);
935:   ksp->data                              = (void*)lgmres;
936:   ksp->ops->buildsolution                = KSPBuildSolution_LGMRES;

938:   ksp->ops->setup                        = KSPSetUp_LGMRES;
939:   ksp->ops->solve                        = KSPSolve_LGMRES;
940:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
941:   ksp->ops->view                         = KSPView_LGMRES;
942:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
943:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
944:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

946:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
947:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
948:                                      KSPGMRESSetPreAllocateVectors_GMRES);
949:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
950:                                     "KSPGMRESSetOrthogonalization_GMRES",
951:                                      KSPGMRESSetOrthogonalization_GMRES);
952:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
953:                                     "KSPGMRESSetRestart_GMRES",
954:                                      KSPGMRESSetRestart_GMRES);
955:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
956:                                     "KSPGMRESSetHapTol_GMRES",
957:                                      KSPGMRESSetHapTol_GMRES);
958:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
959:                                     "KSPGMRESSetCGSRefinementType_GMRES",
960:                                      KSPGMRESSetCGSRefinementType_GMRES);

962:   /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
963:    PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
964:                                      "KSPLGMRESSetConstant_LGMRES",
965:                                       KSPLGMRESSetConstant_LGMRES);

967:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
968:                                     "KSPLGMRESSetAugDim_LGMRES",
969:                                      KSPLGMRESSetAugDim_LGMRES);
970: 

972:   /*defaults */
973:   lgmres->haptol              = 1.0e-30;
974:   lgmres->q_preallocate       = 0;
975:   lgmres->delta_allocate      = LGMRES_DELTA_DIRECTIONS;
976:   lgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
977:   lgmres->nrs                 = 0;
978:   lgmres->sol_temp            = 0;
979:   lgmres->max_k               = LGMRES_DEFAULT_MAXK;
980:   lgmres->Rsvd                = 0;
981:   lgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
982:   lgmres->orthogwork          = 0;
983:   /*LGMRES_MOD - new defaults */
984:   lgmres->aug_dim             = LGMRES_DEFAULT_AUGDIM;
985:   lgmres->aug_ct              = 0; /* start with no aug vectors */
986:   lgmres->approx_constant     = 0;
987:   lgmres->matvecs             = 0;

989:   return(0);
990: }