Actual source code: lgmres.c


  2: #include <../src/ksp/ksp/impls/gmres/lgmres/lgmresimpl.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 PetscErrorCode KSPLGMRESGetNewVectors(KSP, PetscInt);
  8: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP, PetscInt, PetscBool, PetscReal *);
  9: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);

 11: PetscErrorCode KSPLGMRESSetAugDim(KSP ksp, PetscInt dim)
 12: {
 13:   PetscTryMethod((ksp), "KSPLGMRESSetAugDim_C", (KSP, PetscInt), (ksp, dim));
 14:   return 0;
 15: }

 17: PetscErrorCode KSPLGMRESSetConstant(KSP ksp)
 18: {
 19:   PetscTryMethod((ksp), "KSPLGMRESSetConstant_C", (KSP), (ksp));
 20:   return 0;
 21: }

 23: /*
 24:     KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.

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

 29: */
 30: PetscErrorCode KSPSetUp_LGMRES(KSP ksp)
 31: {
 32:   PetscInt    max_k, k, aug_dim;
 33:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

 35:   max_k   = lgmres->max_k;
 36:   aug_dim = lgmres->aug_dim;
 37:   KSPSetUp_GMRES(ksp);

 39:   /* need array of pointers to augvecs*/
 40:   PetscMalloc1(2 * aug_dim + AUG_OFFSET, &lgmres->augvecs);

 42:   lgmres->aug_vecs_allocated = 2 * aug_dim + AUG_OFFSET;

 44:   PetscMalloc1(2 * aug_dim + AUG_OFFSET, &lgmres->augvecs_user_work);
 45:   PetscMalloc1(aug_dim, &lgmres->aug_order);

 47:   /*  for now we will preallocate the augvecs - because aug_dim << restart
 48:      ... also keep in mind that we need to keep augvecs from cycle to cycle*/
 49:   lgmres->aug_vv_allocated = 2 * aug_dim + AUG_OFFSET;
 50:   lgmres->augwork_alloc    = 2 * aug_dim + AUG_OFFSET;

 52:   KSPCreateVecs(ksp, lgmres->aug_vv_allocated, &lgmres->augvecs_user_work[0], 0, NULL);
 53:   PetscMalloc1(max_k + 1, &lgmres->hwork);
 54:   for (k = 0; k < lgmres->aug_vv_allocated; k++) lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
 55:   return 0;
 56: }

 58: /*

 60:     KSPLGMRESCycle - Run lgmres, possibly with restart.  Return residual
 61:                   history if requested.

 63:     input parameters:
 64: .        lgmres  - structure containing parameters and work areas

 66:     output parameters:
 67: .        nres    - residuals (from preconditioned system) at each step.
 68:                   If restarting, consider passing nres+it.  If null,
 69:                   ignored
 70: .        itcount - number of iterations used.   nres[0] to nres[itcount]
 71:                   are defined.  If null, ignored.  If null, ignored.
 72: .        converged - 0 if not converged

 74:     Notes:
 75:     On entry, the value in vector VEC_VV(0) should be
 76:     the initial residual.

 78:  */
 79: PetscErrorCode KSPLGMRESCycle(PetscInt *itcount, KSP ksp)
 80: {
 81:   KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
 82:   PetscReal   res_norm, res;
 83:   PetscReal   hapbnd, tt;
 84:   PetscScalar tmp;
 85:   PetscBool   hapend = PETSC_FALSE;   /* indicates happy breakdown ending */
 86:   PetscInt    loc_it;                 /* local count of # of dir. in Krylov space */
 87:   PetscInt    max_k  = lgmres->max_k; /* max approx space size */
 88:   PetscInt    max_it = ksp->max_it;   /* max # of overall iterations for the method */

 90:   /* LGMRES_MOD - new variables*/
 91:   PetscInt     aug_dim = lgmres->aug_dim;
 92:   PetscInt     spot    = 0;
 93:   PetscInt     order   = 0;
 94:   PetscInt     it_arnoldi; /* number of arnoldi steps to take */
 95:   PetscInt     it_total;   /* total number of its to take (=approx space size)*/
 96:   PetscInt     ii, jj;
 97:   PetscReal    tmp_norm;
 98:   PetscScalar  inv_tmp_norm;
 99:   PetscScalar *avec;

101:   /* Number of pseudo iterations since last restart is the number
102:      of prestart directions */
103:   loc_it = 0;

105:   /* LGMRES_MOD: determine number of arnoldi steps to take */
106:   /* if approx_constant then we keep the space the same size even if
107:      we don't have the full number of aug vectors yet*/
108:   if (lgmres->approx_constant) it_arnoldi = max_k - lgmres->aug_ct;
109:   else it_arnoldi = max_k - aug_dim;

111:   it_total = it_arnoldi + lgmres->aug_ct;

113:   /* initial residual is in VEC_VV(0)  - compute its norm*/
114:   VecNorm(VEC_VV(0), NORM_2, &res_norm);
115:   KSPCheckNorm(ksp, res_norm);
116:   res = res_norm;

118:   /* first entry in right-hand-side of hessenberg system is just
119:      the initial residual norm */
120:   *GRS(0) = res_norm;

122:   /* check for the convergence */
123:   if (!res) {
124:     if (itcount) *itcount = 0;
125:     ksp->reason = KSP_CONVERGED_ATOL;
126:     PetscInfo(ksp, "Converged due to zero residual norm on entry\n");
127:     return 0;
128:   }

130:   /* scale VEC_VV (the initial residual) */
131:   tmp = 1.0 / res_norm;
132:   VecScale(VEC_VV(0), tmp);

134:   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
135:   else ksp->rnorm = 0.0;

137:   /* note: (lgmres->it) is always set one less than (loc_it) It is used in
138:      KSPBUILDSolution_LGMRES, where it is passed to KSPLGMRESBuildSoln.
139:      Note that when KSPLGMRESBuildSoln is called from this function,
140:      (loc_it -1) is passed, so the two are equivalent */
141:   lgmres->it = (loc_it - 1);

143:   /* MAIN ITERATION LOOP BEGINNING*/

145:   /* keep iterating until we have converged OR generated the max number
146:      of directions OR reached the max number of iterations for the method */
147:   (*ksp->converged)(ksp, ksp->its, res, &ksp->reason, ksp->cnvP);

149:   while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
150:     KSPLogResidualHistory(ksp, res);
151:     lgmres->it = (loc_it - 1);
152:     KSPMonitor(ksp, ksp->its, res);

154:     /* see if more space is needed for work vectors */
155:     if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
156:       KSPLGMRESGetNewVectors(ksp, loc_it + 1);
157:       /* (loc_it+1) is passed in as number of the first vector that should
158:           be allocated */
159:     }

161:     /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
162:     if (loc_it < it_arnoldi) { /* Arnoldi */
163:       KSP_PCApplyBAorAB(ksp, VEC_VV(loc_it), VEC_VV(1 + loc_it), VEC_TEMP_MATOP);
164:     } else {                           /*aug step */
165:       order = loc_it - it_arnoldi + 1; /* which aug step */
166:       for (ii = 0; ii < aug_dim; ii++) {
167:         if (lgmres->aug_order[ii] == order) {
168:           spot = ii;
169:           break; /* must have this because there will be duplicates before aug_ct = aug_dim */
170:         }
171:       }

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

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

182:     /* new entry in hessenburg is the 2-norm of our new direction */
183:     VecNorm(VEC_VV(loc_it + 1), NORM_2, &tt);

185:     *HH(loc_it + 1, loc_it)  = tt;
186:     *HES(loc_it + 1, loc_it) = tt;

188:     /* check for the happy breakdown */
189:     hapbnd = PetscAbsScalar(tt / *GRS(loc_it)); /* GRS(loc_it) contains the res_norm from the last iteration  */
190:     if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
191:     if (tt > hapbnd) {
192:       tmp = 1.0 / tt;
193:       VecScale(VEC_VV(loc_it + 1), tmp); /* scale new direction by its norm */
194:     } else {
195:       PetscInfo(ksp, "Detected happy breakdown, current hapbnd = %g tt = %g\n", (double)hapbnd, (double)tt);
196:       hapend = PETSC_TRUE;
197:     }

199:     /* Now apply rotations to new col of hessenberg (and right side of system),
200:        calculate new rotation, and get new residual norm at the same time*/
201:     KSPLGMRESUpdateHessenberg(ksp, loc_it, hapend, &res);
202:     if (ksp->reason) break;

204:     loc_it++;
205:     lgmres->it = (loc_it - 1); /* Add this here in case it has converged */

207:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
208:     ksp->its++;
209:     if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
210:     else ksp->rnorm = 0.0;
211:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

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

215:     /* Catch error in happy breakdown and signal convergence and break from loop */
216:     if (hapend) {
217:       if (!ksp->reason) {
219:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
220:         break;
221:       }
222:     }
223:   }
224:   /* END OF ITERATION LOOP */
225:   KSPLogResidualHistory(ksp, res);

227:   /* Monitor if we know that we will not return for a restart */
228:   if (ksp->reason || ksp->its >= max_it) KSPMonitor(ksp, ksp->its, res);

230:   if (itcount) *itcount = loc_it;

232:   /*
233:     Down here we have to solve for the "best" coefficients of the Krylov
234:     columns, add the solution values together, and possibly unwind the
235:     preconditioning from the solution
236:    */

238:   /* Form the solution (or the solution so far) */
239:   /* Note: must pass in (loc_it-1) for iteration count so that KSPLGMRESBuildSoln
240:      properly navigates */

242:   KSPLGMRESBuildSoln(GRS(0), ksp->vec_sol, ksp->vec_sol, ksp, loc_it - 1);

244:   /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
245:      only if we will be restarting (i.e. this cycle performed it_total
246:      iterations)  */
247:   if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
248:     /*AUG_TEMP contains the new augmentation vector (assigned in  KSPLGMRESBuildSoln) */
249:     if (!lgmres->aug_ct) {
250:       spot = 0;
251:       lgmres->aug_ct++;
252:     } else if (lgmres->aug_ct < aug_dim) {
253:       spot = lgmres->aug_ct;
254:       lgmres->aug_ct++;
255:     } else { /* truncate */
256:       for (ii = 0; ii < aug_dim; ii++) {
257:         if (lgmres->aug_order[ii] == aug_dim) spot = ii;
258:       }
259:     }

261:     VecCopy(AUG_TEMP, AUGVEC(spot));
262:     /*need to normalize */
263:     VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);

265:     inv_tmp_norm = 1.0 / tmp_norm;

267:     VecScale(AUGVEC(spot), inv_tmp_norm);

269:     /*set new aug vector to order 1  - move all others back one */
270:     for (ii = 0; ii < aug_dim; ii++) AUG_ORDER(ii)++;
271:     AUG_ORDER(spot) = 1;

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

276:     /* first do H+*y */
277:     avec = lgmres->hwork;
278:     PetscArrayzero(avec, it_total + 1);
279:     for (ii = 0; ii < it_total + 1; ii++) {
280:       for (jj = 0; jj <= ii + 1 && jj < it_total + 1; jj++) avec[jj] += *HES(jj, ii) * *GRS(ii);
281:     }

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

287:     /*copy answer to aug location  and scale*/
288:     VecCopy(VEC_TEMP, A_AUGVEC(spot));
289:     VecScale(A_AUGVEC(spot), inv_tmp_norm);
290:   }
291:   return 0;
292: }

294: /*
295:     KSPSolve_LGMRES - This routine applies the LGMRES method.

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

300:    Output Parameter:
301: .     outits - number of iterations used

303: */

305: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
306: {
307:   PetscInt    cycle_its; /* iterations done in a call to KSPLGMRESCycle */
308:   PetscInt    itcount;   /* running total of iterations, incl. those in restarts */
309:   KSP_LGMRES *lgmres     = (KSP_LGMRES *)ksp->data;
310:   PetscBool   guess_zero = ksp->guess_zero;
311:   PetscInt    ii; /*LGMRES_MOD variable */


315:   PetscObjectSAWsTakeAccess((PetscObject)ksp);

317:   ksp->its        = 0;
318:   lgmres->aug_ct  = 0;
319:   lgmres->matvecs = 0;

321:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

323:   /* initialize */
324:   itcount = 0;
325:   /*LGMRES_MOD*/
326:   for (ii = 0; ii < lgmres->aug_dim; ii++) lgmres->aug_order[ii] = 0;

328:   while (!ksp->reason) {
329:     /* calc residual - puts in VEC_VV(0) */
330:     KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs);
331:     KSPLGMRESCycle(&cycle_its, ksp);
332:     itcount += cycle_its;
333:     if (itcount >= ksp->max_it) {
334:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
335:       break;
336:     }
337:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
338:   }
339:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
340:   return 0;
341: }

343: /*

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

347: */
348: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
349: {
350:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

352:   PetscFree(lgmres->augvecs);
353:   if (lgmres->augwork_alloc) VecDestroyVecs(lgmres->augwork_alloc, &lgmres->augvecs_user_work[0]);
354:   PetscFree(lgmres->augvecs_user_work);
355:   PetscFree(lgmres->aug_order);
356:   PetscFree(lgmres->hwork);
357:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetConstant_C", NULL);
358:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetAugDim_C", NULL);
359:   KSPDestroy_GMRES(ksp);
360:   return 0;
361: }

363: /*
364:     KSPLGMRESBuildSoln - create the solution from the starting vector and the
365:                       current iterates.

367:     Input parameters:
368:         nrs - work area of size it + 1.
369:         vguess  - index of initial guess
370:         vdest - index of result.  Note that vguess may == vdest (replace
371:                 guess with the solution).
372:         it - HH upper triangular part is a block of size (it+1) x (it+1)

374:      This is an internal routine that knows about the LGMRES internals.
375:  */
376: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar *nrs, Vec vguess, Vec vdest, KSP ksp, PetscInt it)
377: {
378:   PetscScalar tt;
379:   PetscInt    ii, k, j;
380:   KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
381:   /*LGMRES_MOD */
382:   PetscInt it_arnoldi, it_aug;
383:   PetscInt jj, spot = 0;

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

387:   /* If it is < 0, no lgmres steps have been performed */
388:   if (it < 0) {
389:     VecCopy(vguess, vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
390:     return 0;
391:   }

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

395:   /* LGMRES_MOD - determine if we need to use augvecs for the soln  - do not assume that
396:      this is called after the total its allowed for an approx space */
397:   if (lgmres->approx_constant) {
398:     it_arnoldi = lgmres->max_k - lgmres->aug_ct;
399:   } else {
400:     it_arnoldi = lgmres->max_k - lgmres->aug_dim;
401:   }
402:   if (it_arnoldi >= it + 1) {
403:     it_aug     = 0;
404:     it_arnoldi = it + 1;
405:   } else {
406:     it_aug = (it + 1) - it_arnoldi;
407:   }

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

412:   /* solve the upper triangular system - GRS is the right side and HH is
413:      the upper triangular matrix  - put soln in nrs */
415:   if (*HH(it, it) != 0.0) {
416:     nrs[it] = *GRS(it) / *HH(it, it);
417:   } else {
418:     nrs[it] = 0.0;
419:   }

421:   for (ii = 1; ii <= it; ii++) {
422:     k  = it - ii;
423:     tt = *GRS(k);
424:     for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
425:     nrs[k] = tt / *HH(k, k);
426:   }

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

431:   /*LGMRES_MOD - if augmenting has happened we need to form the solution
432:     using the augvecs */
433:   if (!it_aug) { /* all its are from arnoldi */
434:     VecMAXPY(VEC_TEMP, it + 1, nrs, &VEC_VV(0));
435:   } else { /*use aug vecs */
436:     /*first do regular krylov directions */
437:     VecMAXPY(VEC_TEMP, it_arnoldi, nrs, &VEC_VV(0));
438:     /*now add augmented portions - add contribution of aug vectors one at a time*/

440:     for (ii = 0; ii < it_aug; ii++) {
441:       for (jj = 0; jj < lgmres->aug_dim; jj++) {
442:         if (lgmres->aug_order[jj] == (ii + 1)) {
443:           spot = jj;
444:           break; /* must have this because there will be duplicates before aug_ct = aug_dim */
445:         }
446:       }
447:       VecAXPY(VEC_TEMP, nrs[it_arnoldi + ii], AUGVEC(spot));
448:     }
449:   }
450:   /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
451:      preconditioner is "unwound" from right-precondtioning*/
452:   VecCopy(VEC_TEMP, AUG_TEMP);

454:   KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP);

456:   /* add solution to previous solution */
457:   /* put updated solution into vdest.*/
458:   VecCopy(vguess, vdest);
459:   VecAXPY(vdest, 1.0, VEC_TEMP);
460:   return 0;
461: }

463: /*

465:     KSPLGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
466:                             Return new residual.

468:     input parameters:

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

475:     output parameters:
476: .        res - the new residual

478:  */
479: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
480: {
481:   PetscScalar *hh, *cc, *ss, tt;
482:   PetscInt     j;
483:   KSP_LGMRES  *lgmres = (KSP_LGMRES *)(ksp->data);

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

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

494:   for (j = 1; j <= it; j++) {
495:     tt  = *hh;
496:     *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
497:     hh++;
498:     *hh = *cc++ * *hh - (*ss++ * tt);
499:     /* hh, cc, and ss have all been incremented one by end of loop */
500:   }

502:   /*
503:     compute the new plane rotation, and apply it to:
504:      1) the right-hand-side of the Hessenberg system (GRS)
505:         note: it affects GRS(it) and GRS(it+1)
506:      2) the new column of the Hessenberg matrix
507:         note: it affects HH(it,it) which is currently pointed to
508:         by hh and HH(it+1, it) (*(hh+1))
509:     thus obtaining the updated value of the residual...
510:   */

512:   /* compute new plane rotation */

514:   if (!hapend) {
515:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
516:     if (tt == 0.0) {
517:       ksp->reason = KSP_DIVERGED_NULL;
518:       return 0;
519:     }
520:     *cc = *hh / tt;       /* new cosine value */
521:     *ss = *(hh + 1) / tt; /* new sine value */

523:     /* apply to 1) and 2) */
524:     *GRS(it + 1) = -(*ss * *GRS(it));
525:     *GRS(it)     = PetscConj(*cc) * *GRS(it);
526:     *hh          = PetscConj(*cc) * *hh + *ss * *(hh + 1);

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

531:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
532:             another rotation matrix (so RH doesn't change).  The new residual is
533:             always the new sine term times the residual from last time (GRS(it)),
534:             but now the new sine rotation would be zero...so the residual should
535:             be zero...so we will multiply "zero" by the last residual.  This might
536:             not be exactly what we want to do here -could just return "zero". */

538:     *res = 0.0;
539:   }
540:   return 0;
541: }

543: /*

545:    KSPLGMRESGetNewVectors - This routine allocates more work vectors, starting from
546:                          VEC_VV(it)

548: */
549: static PetscErrorCode KSPLGMRESGetNewVectors(KSP ksp, PetscInt it)
550: {
551:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
552:   PetscInt    nwork  = lgmres->nwork_alloc; /* number of work vector chunks allocated */
553:   PetscInt    nalloc;                       /* number to allocate */
554:   PetscInt    k;

556:   nalloc = lgmres->delta_allocate; /* number of vectors to allocate
557:                                       in a single chunk */

559:   /* Adjust the number to allocate to make sure that we don't exceed the
560:      number of available slots (lgmres->vecs_allocated)*/
561:   if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated) nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
562:   if (!nalloc) return 0;

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

566:   /* work vectors */
567:   KSPCreateVecs(ksp, nalloc, &lgmres->user_work[nwork], 0, NULL);
568:   /* specify size of chunk allocated */
569:   lgmres->mwork_alloc[nwork] = nalloc;

571:   for (k = 0; k < nalloc; k++) lgmres->vecs[it + VEC_OFFSET + k] = lgmres->user_work[nwork][k];

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

575:   /* increment the number of work vector chunks */
576:   lgmres->nwork_alloc++;
577:   return 0;
578: }

580: /*

582:    KSPBuildSolution_LGMRES

584:      Input Parameter:
585: .     ksp - the Krylov space object
586: .     ptr-

588:    Output Parameter:
589: .     result - the solution

591:    Note: this calls KSPLGMRESBuildSoln - the same function that KSPLGMRESCycle
592:    calls directly.

594: */
595: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp, Vec ptr, Vec *result)
596: {
597:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

599:   if (!ptr) {
600:     if (!lgmres->sol_temp) { VecDuplicate(ksp->vec_sol, &lgmres->sol_temp); }
601:     ptr = lgmres->sol_temp;
602:   }
603:   if (!lgmres->nrs) {
604:     /* allocate the work area */
605:     PetscMalloc1(lgmres->max_k, &lgmres->nrs);
606:   }

608:   KSPLGMRESBuildSoln(lgmres->nrs, ksp->vec_sol, ptr, ksp, lgmres->it);
609:   if (result) *result = ptr;
610:   return 0;
611: }

613: PetscErrorCode KSPView_LGMRES(KSP ksp, PetscViewer viewer)
614: {
615:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
616:   PetscBool   iascii;

618:   KSPView_GMRES(ksp, viewer);
619:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
620:   if (iascii) {
621:     /*LGMRES_MOD */
622:     PetscViewerASCIIPrintf(viewer, "  aug. dimension=%" PetscInt_FMT "\n", lgmres->aug_dim);
623:     if (lgmres->approx_constant) PetscViewerASCIIPrintf(viewer, "  approx. space size was kept constant.\n");
624:     PetscViewerASCIIPrintf(viewer, "  number of matvecs=%" PetscInt_FMT "\n", lgmres->matvecs);
625:   }
626:   return 0;
627: }

629: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
630: {
631:   PetscInt    aug;
632:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
633:   PetscBool   flg    = PETSC_FALSE;

635:   KSPSetFromOptions_GMRES(ksp, PetscOptionsObject);
636:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP LGMRES Options");
637:   PetscOptionsBool("-ksp_lgmres_constant", "Use constant approx. space size", "KSPGMRESSetConstant", lgmres->approx_constant, &lgmres->approx_constant, NULL);
638:   PetscOptionsInt("-ksp_lgmres_augment", "Number of error approximations to augment the Krylov space with", "KSPLGMRESSetAugDim", lgmres->aug_dim, &aug, &flg);
639:   if (flg) KSPLGMRESSetAugDim(ksp, aug);
640:   PetscOptionsHeadEnd();
641:   return 0;
642: }

644: /*functions for extra lgmres options here*/
645: static PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)
646: {
647:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

649:   lgmres->approx_constant = PETSC_TRUE;
650:   return 0;
651: }

653: static PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp, PetscInt aug_dim)
654: {
655:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

659:   lgmres->aug_dim = aug_dim;
660:   return 0;
661: }

663: /* end new lgmres functions */

665: /*MC
666:     KSPLGMRES - Augments the standard GMRES approximation space with approximations to
667:                 the error from previous restart cycles.

669:   Options Database Keys:
670: +   -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
671: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
672: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
673:                             vectors are allocated as needed)
674: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
675: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
676: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
677:                                   stability of the classical Gram-Schmidt  orthogonalization.
678: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
679: .   -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
680: -   -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)

682:     To run LGMRES(m, k) as described in the above paper, use:
683:        -ksp_gmres_restart <m+k>
684:        -ksp_lgmres_augment <k>

686:   Level: beginner

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

691:    References:
692: .  * - A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for accelerating the convergence of restarted GMRES. SIAM Journal on Matrix Analysis and Applications, 26 (2005).

694:   Developer Notes:
695:     This object is subclassed off of KSPGMRES

697:   Contributed by: Allison Baker

699: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPGMRES`,
700:           `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
701:           `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
702:           `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPLGMRESSetAugDim()`,
703:           `KSPGMRESSetConstant()`

705: M*/

707: PETSC_EXTERN PetscErrorCode KSPCreate_LGMRES(KSP ksp)
708: {
709:   KSP_LGMRES *lgmres;

711:   PetscNew(&lgmres);

713:   ksp->data               = (void *)lgmres;
714:   ksp->ops->buildsolution = KSPBuildSolution_LGMRES;

716:   ksp->ops->setup                        = KSPSetUp_LGMRES;
717:   ksp->ops->solve                        = KSPSolve_LGMRES;
718:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
719:   ksp->ops->view                         = KSPView_LGMRES;
720:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
721:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
722:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

724:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
725:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
726:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

728:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES);
729:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES);
730:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", KSPGMRESGetOrthogonalization_GMRES);
731:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES);
732:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_GMRES);
733:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES);
734:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES);
735:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", KSPGMRESGetCGSRefinementType_GMRES);

737:   /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
738:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetConstant_C", KSPLGMRESSetConstant_LGMRES);
739:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetAugDim_C", KSPLGMRESSetAugDim_LGMRES);

741:   /*defaults */
742:   lgmres->haptol         = 1.0e-30;
743:   lgmres->q_preallocate  = 0;
744:   lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
745:   lgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
746:   lgmres->nrs            = NULL;
747:   lgmres->sol_temp       = NULL;
748:   lgmres->max_k          = LGMRES_DEFAULT_MAXK;
749:   lgmres->Rsvd           = NULL;
750:   lgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
751:   lgmres->orthogwork     = NULL;

753:   /*LGMRES_MOD - new defaults */
754:   lgmres->aug_dim         = LGMRES_DEFAULT_AUGDIM;
755:   lgmres->aug_ct          = 0; /* start with no aug vectors */
756:   lgmres->approx_constant = PETSC_FALSE;
757:   lgmres->matvecs         = 0;
758:   return 0;
759: }