Actual source code: lgmres.c

 2:  #include lgmresp.h

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

 13: int KSPLGMRESSetAugDim(KSP ksp, int dim)
 14: {

 18:   PetscTryMethod((ksp),KSPLGMRESSetAugDim_C,(KSP,int),(ksp,dim));
 19:   return(0);
 20: }

 24: int KSPLGMRESSetConstant(KSP ksp)
 25: {

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

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

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

 39: */
 42: int    KSPSetUp_LGMRES(KSP ksp)
 43: {
 44:   unsigned  int size,hh,hes,rs,cc;
 45:   int           ierr,max_k,k, aug_dim;
 46:   KSP_LGMRES    *lgmres = (KSP_LGMRES *)ksp->data;

 49:   if (ksp->pc_side == PC_SYMMETRIC) {
 50:     SETERRQ(2,"no symmetric preconditioning for KSPLGMRES");
 51:   }



 55:   max_k         = lgmres->max_k;
 56:   aug_dim       = lgmres->aug_dim;
 57:   hh            = (max_k + 2) * (max_k + 1);
 58:   hes           = (max_k + 1) * (max_k + 1);
 59:   rs            = (max_k + 2);
 60:   cc            = (max_k + 1);  /* SS and CC are the same size */
 61:   size          = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);

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

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

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

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

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

130:   return(0);
131: }


134: /*

136:     LGMRESCycle - Run lgmres, possibly with restart.  Return residual 
137:                   history if requested.

139:     input parameters:
140: .         lgmres  - structure containing parameters and work areas

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

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


156:  */
159: int LGMREScycle(int *itcount,KSP ksp)
160: {

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


185:   /* Number of pseudo iterations since last restart is the number 
186:      of prestart directions */
187:   loc_it = 0;

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

198:   it_total =  it_arnoldi + lgmres->aug_ct;

200:   /* initial residual is in VEC_VV(0)  - compute its norm*/
201:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
202:   res    = res_norm;
203: 
204:   /* first entry in right-hand-side of hessenberg system is just 
205:      the initial residual norm */
206:   *GRS(0) = res_norm;

208:  /* check for the convergence */
209:   if (!res) {
210:      if (itcount) *itcount = 0;
211:      ksp->reason = KSP_CONVERGED_ATOL;
212:      PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entry\n");
213:      return(0);
214:   }

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

219:   /* FYI: AMS calls are for memory snooper */
220:   PetscObjectTakeAccess(ksp);
221:   ksp->rnorm = res;
222:   PetscObjectGrantAccess(ksp);


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

231: 
232:   /* MAIN ITERATION LOOP BEGINNING*/


235:   /* keep iterating until we have converged OR generated the max number
236:      of directions OR reached the max number of iterations for the method */
237:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
238: 
239:   while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
240:      KSPLogResidualHistory(ksp,res);
241:      lgmres->it = (loc_it - 1);
242:      KSPMonitor(ksp,ksp->its,res);

244:     /* see if more space is needed for work vectors */
245:     if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
246:        LGMRESGetNewVectors(ksp,loc_it+1);
247:       /* (loc_it+1) is passed in as number of the first vector that should
248:          be allocated */
249:     }

251:     /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
252:     if (loc_it < it_arnoldi) { /* arnoldi */
253:        KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
254:     } else { /*aug step */
255:        order = loc_it - it_arnoldi + 1; /* which aug step */
256:        for (ii=0; ii<aug_dim; ii++) {
257:            if (lgmres->aug_order[ii] == order) {
258:               spot = ii;
259:               break; /* must have this because there will be duplicates before aug_ct = aug_dim */
260:             }
261:         }

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

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

272:     /* new entry in hessenburg is the 2-norm of our new direction */
273:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
274:     *HH(loc_it+1,loc_it)   = tt;
275:     *HES(loc_it+1,loc_it)  = tt;


278:     /* check for the happy breakdown */
279:     hapbnd  = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration  */
280:     if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
281:     if (tt > hapbnd) {
282:        tmp = 1.0/tt;
283:        VecScale(&tmp,VEC_VV(loc_it+1)); /* scale new direction by its norm */
284:     } else {
285:        PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",hapbnd,tt);
286:        hapend = PETSC_TRUE;
287:     }

289:     /* Now apply rotations to new col of hessenberg (and right side of system), 
290:        calculate new rotation, and get new residual norm at the same time*/
291:     LGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
292:     loc_it++;
293:     lgmres->it  = (loc_it-1);  /* Add this here in case it has converged */
294: 
295:     PetscObjectTakeAccess(ksp);
296:     ksp->its++;
297:     ksp->rnorm = res;
298:     PetscObjectGrantAccess(ksp);

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

302:     /* Catch error in happy breakdown and signal convergence and break from loop */
303:     if (hapend) {
304:       if (!ksp->reason) {
305:         SETERRQ1(0,"You reached the happy break down,but convergence was not indicated. Residual norm = %g",res);
306:       }
307:       break;
308:     }
309:   }
310:   /* END OF ITERATION LOOP */

312:   KSPLogResidualHistory(ksp,res);

314:   /* Monitor if we know that we will not return for a restart */
315:   if (ksp->reason || ksp->its >= max_it) {
316:     KSPMonitor(ksp, ksp->its, res);
317:   }

319:   if (itcount) *itcount    = loc_it;

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

331:   BuildLgmresSoln(GRS(0),VEC_SOLN,VEC_SOLN,ksp,loc_it-1);


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

339:      /*AUG_TEMP contains the new augmentation vector (assigned in  BuildLgmresSoln) */
340:     if (lgmres->aug_ct == 0) {
341:         spot = 0;
342:         lgmres->aug_ct++;
343:      } else if (lgmres->aug_ct < aug_dim) {
344:         spot = lgmres->aug_ct;
345:         lgmres->aug_ct++;
346:      } else { /* truncate */
347:         for (ii=0; ii<aug_dim; ii++) {
348:            if (lgmres->aug_order[ii] == aug_dim) {
349:               spot = ii;
350:             }
351:         }
352:      }

354: 

356:      VecCopy(AUG_TEMP, AUGVEC(spot));
357:      /*need to normalize */
358:      VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
359:      inv_tmp_norm = 1.0/tmp_norm;
360:      VecScale(&inv_tmp_norm, AUGVEC(spot));

362:      /*set new aug vector to order 1  - move all others back one */
363:      for (ii=0; ii < aug_dim; ii++) {
364:         AUG_ORDER(ii)++;
365:      }
366:      AUG_ORDER(spot) = 1;

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

371: 
372:      /* first do H+*y */
373:      VecSet(&zero, AUG_TEMP);
374:      VecGetArray(AUG_TEMP, &avec);
375:      for (ii=0; ii < it_total + 1; ii++) {
376:         for (jj=0; jj <= ii+1; jj++) {
377:            avec[jj] += *HES(jj ,ii) * *GRS(ii);
378:         }
379:      }

381:      /*now multiply result by V+ */
382:      VecSet(&zero, VEC_TEMP);
383:      VecMAXPY(it_total+1, avec, VEC_TEMP, &VEC_VV(0)); /*answer is in VEC_TEMP*/
384:      VecRestoreArray(AUG_TEMP, &avec);
385: 
386:      /*copy answer to aug location  and scale*/
387:      VecCopy(VEC_TEMP,  A_AUGVEC(spot));
388:      VecScale(&inv_tmp_norm, A_AUGVEC(spot));


391:   }
392:   return(0);
393: }

395: /*  
396:     KSPSolve_LGMRES - This routine applies the LGMRES method.


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

402:    Output Parameter:
403: .     outits - number of iterations used

405: */

409: int KSPSolve_LGMRES(KSP ksp)
410: {
411:   int        ierr;
412:   int        cycle_its; /* iterations done in a call to LGMREScycle */
413:   int        itcount;   /* running total of iterations, incl. those in restarts */
414:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
415:   PetscTruth guess_zero = ksp->guess_zero;
416:   int        ii;        /*LGMRES_MOD variable */

419:   if (ksp->calc_sings && !lgmres->Rsvd) {
420:      SETERRQ(1,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
421:   }
422:   PetscObjectTakeAccess(ksp);
423:   ksp->its        = 0;
424:   lgmres->aug_ct  = 0;
425:   lgmres->matvecs = 0;
426:   PetscObjectGrantAccess(ksp);

428:   /* initialize */
429:   itcount  = 0;
430:   ksp->reason = KSP_CONVERGED_ITERATING;
431:   /*LGMRES_MOD*/
432:   for (ii=0; ii<lgmres->aug_dim; ii++) {
433:      lgmres->aug_order[ii] = 0;
434:   }

436:   while (!ksp->reason) {
437:      /* calc residual - puts in VEC_VV(0) */
438:     KSPInitialResidual(ksp,VEC_SOLN,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),VEC_RHS);
439:     LGMREScycle(&cycle_its,ksp);
440:     itcount += cycle_its;
441:     if (itcount >= ksp->max_it) {
442:       ksp->reason = KSP_DIVERGED_ITS;
443:       break;
444:     }
445:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
446:   }
447:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
448:   return(0);
449: }

451: /*

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

455: */
458: int KSPDestroy_LGMRES(KSP ksp)
459: {
460:   KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
461:   int        i,ierr;

464:   /* Free the Hessenberg matrices */
465:   if (lgmres->hh_origin) {PetscFree(lgmres->hh_origin);}

467:   /* Free pointers to user variables */
468:   if (lgmres->vecs) {PetscFree(lgmres->vecs);}

470:   /*LGMRES_MOD - free pointers for extra vectors */
471:   if (lgmres->augvecs) {PetscFree(lgmres->augvecs);}

473:   /* free work vectors */
474:   for (i=0; i < lgmres->nwork_alloc; i++) {
475:     VecDestroyVecs(lgmres->user_work[i],lgmres->mwork_alloc[i]);
476:   }
477:   if (lgmres->user_work)  {PetscFree(lgmres->user_work);}

479:   /*LGMRES_MOD - free aug work vectors also */
480:   /*this was all allocated as one "chunk" */
481:   VecDestroyVecs(lgmres->augvecs_user_work[0],lgmres->augwork_alloc);
482:   if (lgmres->augvecs_user_work)  {PetscFree(lgmres->augvecs_user_work);}
483:   if (lgmres->aug_order) {PetscFree(lgmres->aug_order);}

485:   if (lgmres->mwork_alloc) {PetscFree(lgmres->mwork_alloc);}
486:   if (lgmres->nrs) {PetscFree(lgmres->nrs);}
487:   if (lgmres->sol_temp) {VecDestroy(lgmres->sol_temp);}
488:   if (lgmres->Rsvd) {PetscFree(lgmres->Rsvd);}
489:   if (lgmres->Dsvd) {PetscFree(lgmres->Dsvd);}
490:   PetscFree(lgmres);
491:   return(0);
492: }

494: /*
495:     BuildLgmresSoln - create the solution from the starting vector and the
496:                       current iterates.

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

505:      This is an internal routine that knows about the LGMRES internals.
506:  */
509: static int BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,int it)
510: {
511:   PetscScalar  tt,zero = 0.0,one = 1.0;
512:   int          ierr,ii,k,j;
513:   KSP_LGMRES   *lgmres = (KSP_LGMRES *)(ksp->data);
514:   /*LGMRES_MOD */
515:   int          it_arnoldi, it_aug;
516:   int          jj, spot = 0;

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

521:   /* If it is < 0, no lgmres steps have been performed */
522:   if (it < 0) {
523:     if (vdest != vguess) {
524:       VecCopy(vguess,vdest);
525:     }
526:     return(0);
527:   }

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

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

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

548: 
549:   /* solve the upper triangular system - GRS is the right side and HH is 
550:      the upper triangular matrix  - put soln in nrs */
551:   if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
552:   if (*HH(it,it) != 0.0) {
553:      nrs[it] = *GRS(it) / *HH(it,it);
554:   } else {
555:      nrs[it] = 0.0;
556:   }

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

565:   /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
566:   VecSet(&zero,VEC_TEMP); /* set VEC_TEMP components to 0 */

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


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

592:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);

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

601:   return(0);
602: }

604: /*

606:     LGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.  
607:                             Return new residual.

609:     input parameters:

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

616:     output parameters:
617: .        res - the new residual
618:         
619:  */
622: static int LGMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
623: {
624:   PetscScalar   *hh,*cc,*ss,tt;
625:   int           j;
626:   KSP_LGMRES    *lgmres = (KSP_LGMRES *)(ksp->data);

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

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

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

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

660:   /* compute new plane rotation */

662:   if (!hapend) {
663: #if defined(PETSC_USE_COMPLEX)
664:     tt        = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
665: #else
666:     tt        = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
667: #endif
668:     if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
669:     *cc       = *hh / tt;   /* new cosine value */
670:     *ss       = *(hh+1) / tt;  /* new sine value */

672:     /* apply to 1) and 2) */
673:     *GRS(it+1) = - (*ss * *GRS(it));
674: #if defined(PETSC_USE_COMPLEX)
675:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
676:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);
677: #else
678:     *GRS(it)   = *cc * *GRS(it);
679:     *hh        = *cc * *hh + *ss * *(hh+1);
680: #endif

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

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

697: /*

699:    LGMRESGetNewVectors - This routine allocates more work vectors, starting from 
700:                          VEC_VV(it) 
701:                          
702: */
705: static int LGMRESGetNewVectors(KSP ksp,int it)
706: {
707:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
708:   int        nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
709:   int        nalloc;                      /* number to allocate */
710:   int        k,ierr;
711: 
713:   nalloc = lgmres->delta_allocate; /* number of vectors to allocate 
714:                                       in a single chunk */

716:   /* Adjust the number to allocate to make sure that we don't exceed the
717:      number of available slots (lgmres->vecs_allocated)*/
718:   if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
719:     nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
720:   }
721:   if (!nalloc) return(0);

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

725:   /* work vectors */
726:   VecDuplicateVecs(ksp->vec_rhs,nalloc,&lgmres->user_work[nwork]);
727:   PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
728:   /* specify size of chunk allocated */
729:   lgmres->mwork_alloc[nwork] = nalloc;

731:   for (k=0; k < nalloc; k++) {
732:     lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
733:   }
734: 

736:   /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
737: 

739:   /* increment the number of work vector chunks */
740:   lgmres->nwork_alloc++;
741:   return(0);
742: }

744: /* 

746:    KSPBuildSolution_LGMRES

748:      Input Parameter:
749: .     ksp - the Krylov space object
750: .     ptr-

752:    Output Parameter:
753: .     result - the solution

755:    Note: this calls BuildLgmresSoln - the same function that LGMREScycle
756:    calls directly.  

758: */
761: int KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
762: {
763:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
764:   int        ierr;

767:   if (!ptr) {
768:     if (!lgmres->sol_temp) {
769:       VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
770:       PetscLogObjectParent(ksp,lgmres->sol_temp);
771:     }
772:     ptr = lgmres->sol_temp;
773:   }
774:   if (!lgmres->nrs) {
775:     /* allocate the work area */
776:     PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
777:     PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
778:   }
779: 
780:   BuildLgmresSoln(lgmres->nrs,VEC_SOLN,ptr,ksp,lgmres->it);
781:   *result = ptr;
782: 
783:   return(0);
784: }

786: /*

788:    KSPView_LGMRES -Prints information about the current Krylov method 
789:                   being used.

791:  */
794: int KSPView_LGMRES(KSP ksp,PetscViewer viewer)
795: {
796:   KSP_LGMRES   *lgmres = (KSP_LGMRES *)ksp->data;
797:   const char   *cstr;
798:   int          ierr;
799:   PetscTruth   isascii,isstring;

802:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
803:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
804:   if (lgmres->orthog == KSPGMRESClassicalGramSchmidtOrthogonalization) {
805:     if (lgmres->cgstype == KSP_GMRES_CGS_REFINE_NEVER) {
806:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement";
807:     } else if (lgmres->cgstype == KSP_GMRES_CGS_REFINE_ALWAYS) {
808:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement";
809:     } else {
810:       cstr = "Classical (unmodified) Gram-Schmidt Orthogonalization with one step of iterative refinement when needed";
811:     }
812:   } else if (lgmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
813:     cstr = "Modified Gram-Schmidt Orthogonalization";
814:   } else {
815:     cstr = "unknown orthogonalization";
816:   }
817:   if (isascii) {
818:     PetscViewerASCIIPrintf(viewer,"  LGMRES: restart=%d, using %s\n",lgmres->max_k,cstr);
819:     /*LGMRES_MOD */
820:     PetscViewerASCIIPrintf(viewer,"  LGMRES: aug. dimension=%d\n",lgmres->aug_dim);
821:     if (lgmres->approx_constant) {
822:        PetscViewerASCIIPrintf(viewer,"  LGMRES: approx. space size was kept constant.\n");
823:     }
824:     PetscViewerASCIIPrintf(viewer,"  LGMRES: number of matvecs=%d\n",lgmres->matvecs);

826:     PetscViewerASCIIPrintf(viewer,"  LGMRES: happy breakdown tolerance %g\n",lgmres->haptol);
827:   } else if (isstring) {
828:     PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,lgmres->max_k);
829:   } else {
830:     SETERRQ1(1,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
831:   }
832:   return(0);


835: }

839: int KSPSetFromOptions_LGMRES(KSP ksp)
840: {
841:   int         ierr, restart, aug,indx;
842:   PetscReal   haptol;
843:   KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
844:   PetscTruth  flg;
845:   const char  *types[] = {"never","ifneeded","always"};

848:   PetscOptionsHead("KSP LGMRES Options");

850: 

852:     PetscOptionsInt("-ksp_gmres_restart","For LGMRES, this is the maximum size of the approximation space","KSPGMRESSetRestart",lgmres->max_k,&restart,&flg);
853:     if (flg) { KSPGMRESSetRestart(ksp,restart); }
854:     PetscOptionsReal("-ksp_gmres_haptol","Tolerance for declaring exact convergence (happy ending)","KSPGMRESSetHapTol",lgmres->haptol,&haptol,&flg);
855:     if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
856:     PetscOptionsName("-ksp_gmres_preallocate","Preallocate all Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
857:     if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
858:     PetscOptionsLogicalGroupBegin("-ksp_gmres_classicalgramschmidt","Use classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
859:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);}
860:     PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Use modified Gram-Schmidt (slow but more stable)","KSPGMRESSetOrthogonalization",&flg);
861:     if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
862:     PetscOptionsEList("-ksp_gmres_cgs_refinement_type","Type of iterative refinement for classical (unmodified) Gram-Schmidt","KSPGMRESSetCGSRefinementType()",types,3,types[lgmres->cgstype],&indx,&flg);
863:     if (flg) {
864:       KSPGMRESSetCGSRefinementType(ksp,(KSPGMRESCGSRefinementType)indx);
865:     }

867:     PetscOptionsName("-ksp_gmres_krylov_monitor","Graphically plot the Krylov directions","KSPSetMonitor",&flg);
868:     if (flg) {
869:       PetscViewers viewers;
870:       PetscViewersCreate(ksp->comm,&viewers);
871:       KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
872:     }

874: /* LGMRES_MOD - specify number of augmented vectors and whether the space should be a constant size*/
875:      PetscOptionsName("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",&flg);
876:     /*if (flg) {KSPGMRESSetConstant(ksp);}*/ /*<--doesn't like this */
877:     if (flg) { lgmres->approx_constant = 1; }                     /* in favor of this line....*/

879:     PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
880:     if (flg) { KSPLGMRESSetAugDim(ksp,aug); }



884:   PetscOptionsTail();
885:   return(0);
886: }


889: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
890: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);

892: /*functions for extra lgmres options here*/
893: EXTERN_C_BEGIN
896: int KSPLGMRESSetConstant_LGMRES(KSP ksp)
897: {
898:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
900:   lgmres->approx_constant = 1;
901: 
902:   return(0);
903: }
904: EXTERN_C_END

906: EXTERN_C_BEGIN
909: int KSPLGMRESSetAugDim_LGMRES(KSP ksp,int aug_dim)
910: {
911:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;


915:   if (aug_dim < 0) SETERRQ(1,"Augmentation dimension must be positive");
916:   if (aug_dim > (lgmres->max_k -1))  SETERRQ(1,"Augmentation dimension must be <= (restart size-1)");

918:   lgmres->aug_dim = aug_dim;

920:   return(0);
921: }
922: EXTERN_C_END


925: /* end new lgmres functions */


928: /* use these options from gmres */
929: EXTERN_C_BEGIN
930: EXTERN int KSPGMRESSetHapTol_GMRES(KSP,double);
931: EXTERN int KSPGMRESSetPreAllocateVectors_GMRES(KSP);
932: EXTERN int KSPGMRESSetRestart_GMRES(KSP,int);
933: EXTERN int KSPGMRESSetOrthogonalization_GMRES(KSP,int (*)(KSP,int));
934: EXTERN int KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);
935: EXTERN_C_END

937: /*MC
938:      KSPLGMRES - Augments the standard GMRES approximation space with approximation to
939:                  the error from previous restart cycles.

941:    Options Database Keys:
942: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
943: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
944: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 
945:                              vectors are allocated as needed)
946: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
947: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
948: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 
949:                                    stability of the classical Gram-Schmidt  orthogonalization.
950: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
951: .   -ksp_lgmres_constant - Use constant approx. space size
952: -   -ksp_lgmres_augment <n> - Number of error approximations to augment the Krylov space with

954:     Described in:
955:      A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
956:      accelerating the convergence of restarted GMRES. Submitted to SIAM
957:      Journal on Matrix Analysis and Applications. Also available as
958:      Technical Report #CU-CS-945-03, University of Colorado, Department of
959:      Computer Science, January, 2003. 

961:    Level: beginner

963:    Contributed by: Allison Baker

965: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
966:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
967:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
968:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor(), KSPLGMRESSetAugDim(),
969:            KSPGMRESSetConstant()

971: M*/

973: EXTERN_C_BEGIN
976: int KSPCreate_LGMRES(KSP ksp)
977: {
978:   KSP_LGMRES *lgmres;
979:   int        ierr;

982:   PetscNew(KSP_LGMRES,&lgmres);
983:   PetscMemzero(lgmres,sizeof(KSP_LGMRES));
984:   PetscLogObjectMemory(ksp,sizeof(KSP_LGMRES));
985:   ksp->data                              = (void*)lgmres;
986:   ksp->ops->buildsolution                = KSPBuildSolution_LGMRES;

988:   ksp->ops->setup                        = KSPSetUp_LGMRES;
989:   ksp->ops->solve                        = KSPSolve_LGMRES;
990:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
991:   ksp->ops->view                         = KSPView_LGMRES;
992:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
993:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
994:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

996:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
997:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
998:                                      KSPGMRESSetPreAllocateVectors_GMRES);
999:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
1000:                                     "KSPGMRESSetOrthogonalization_GMRES",
1001:                                      KSPGMRESSetOrthogonalization_GMRES);
1002:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
1003:                                     "KSPGMRESSetRestart_GMRES",
1004:                                      KSPGMRESSetRestart_GMRES);
1005:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
1006:                                     "KSPGMRESSetHapTol_GMRES",
1007:                                      KSPGMRESSetHapTol_GMRES);
1008:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
1009:                                     "KSPGMRESSetCGSRefinementType_GMRES",
1010:                                      KSPGMRESSetCGSRefinementType_GMRES);

1012:   /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
1013:    PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
1014:                                      "KSPLGMRESSetConstant_LGMRES",
1015:                                       KSPLGMRESSetConstant_LGMRES);

1017:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
1018:                                     "KSPLGMRESSetAugDim_LGMRES",
1019:                                      KSPLGMRESSetAugDim_LGMRES);
1020: 

1022:   /*defaults */
1023:   lgmres->haptol              = 1.0e-30;
1024:   lgmres->q_preallocate       = 0;
1025:   lgmres->delta_allocate      = LGMRES_DELTA_DIRECTIONS;
1026:   lgmres->orthog              = KSPGMRESClassicalGramSchmidtOrthogonalization;
1027:   lgmres->nrs                 = 0;
1028:   lgmres->sol_temp            = 0;
1029:   lgmres->max_k               = LGMRES_DEFAULT_MAXK;
1030:   lgmres->Rsvd                = 0;
1031:   lgmres->cgstype             = KSP_GMRES_CGS_REFINE_NEVER;
1032:   /*LGMRES_MOD - new defaults */
1033:   lgmres->aug_dim             = LGMRES_DEFAULT_AUGDIM;
1034:   lgmres->aug_ct              = 0; /* start with no aug vectors */
1035:   lgmres->approx_constant     = 0;
1036:   lgmres->matvecs             = 0;

1038:   return(0);
1039: }
1040: EXTERN_C_END