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: /*
 59:     KSPLGMRESCycle - Run lgmres, possibly with restart.  Return residual
 60:                   history if requested.

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

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

 73:     Notes:
 74:     On entry, the value in vector VEC_VV(0) should be
 75:     the initial residual.
 76:  */
 77: PetscErrorCode KSPLGMRESCycle(PetscInt *itcount, KSP ksp)
 78: {
 79:   KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
 80:   PetscReal   res_norm, res;
 81:   PetscReal   hapbnd, tt;
 82:   PetscScalar tmp;
 83:   PetscBool   hapend = PETSC_FALSE;   /* indicates happy breakdown ending */
 84:   PetscInt    loc_it;                 /* local count of # of dir. in Krylov space */
 85:   PetscInt    max_k  = lgmres->max_k; /* max approx space size */
 86:   PetscInt    max_it = ksp->max_it;   /* max # of overall iterations for the method */

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

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

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

109:   it_total = it_arnoldi + lgmres->aug_ct;

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

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

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

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

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

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

141:   /* MAIN ITERATION LOOP BEGINNING*/

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

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

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

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

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

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

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

183:     *HH(loc_it + 1, loc_it)  = tt;
184:     *HES(loc_it + 1, loc_it) = tt;

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

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

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

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

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

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

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

228:   if (itcount) *itcount = loc_it;

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

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

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

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

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

263:     inv_tmp_norm = 1.0 / tmp_norm;

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

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

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

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

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

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

292: /*
293:     KSPSolve_LGMRES - This routine applies the LGMRES method.

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

298:    Output Parameter:
299: .     outits - number of iterations used

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


312:   PetscObjectSAWsTakeAccess((PetscObject)ksp);

314:   ksp->its        = 0;
315:   lgmres->aug_ct  = 0;
316:   lgmres->matvecs = 0;

318:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

320:   /* initialize */
321:   itcount = 0;
322:   /*LGMRES_MOD*/
323:   for (ii = 0; ii < lgmres->aug_dim; ii++) lgmres->aug_order[ii] = 0;

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

340: /*
341:    KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
342: */
343: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
344: {
345:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

347:   PetscFree(lgmres->augvecs);
348:   if (lgmres->augwork_alloc) VecDestroyVecs(lgmres->augwork_alloc, &lgmres->augvecs_user_work[0]);
349:   PetscFree(lgmres->augvecs_user_work);
350:   PetscFree(lgmres->aug_order);
351:   PetscFree(lgmres->hwork);
352:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetConstant_C", NULL);
353:   PetscObjectComposeFunction((PetscObject)ksp, "KSPLGMRESSetAugDim_C", NULL);
354:   KSPDestroy_GMRES(ksp);
355:   return 0;
356: }

358: /*
359:     KSPLGMRESBuildSoln - create the solution from the starting vector and the
360:                       current iterates.

362:     Input parameters:
363:         nrs - work area of size it + 1.
364:         vguess  - index of initial guess
365:         vdest - index of result.  Note that vguess may == vdest (replace
366:                 guess with the solution).
367:         it - HH upper triangular part is a block of size (it+1) x (it+1)

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

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

382:   /* If it is < 0, no lgmres steps have been performed */
383:   if (it < 0) {
384:     VecCopy(vguess, vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
385:     return 0;
386:   }

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

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

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

407:   /* solve the upper triangular system - GRS is the right side and HH is
408:      the upper triangular matrix  - put soln in nrs */
410:   if (*HH(it, it) != 0.0) {
411:     nrs[it] = *GRS(it) / *HH(it, it);
412:   } else {
413:     nrs[it] = 0.0;
414:   }

416:   for (ii = 1; ii <= it; ii++) {
417:     k  = it - ii;
418:     tt = *GRS(k);
419:     for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
420:     nrs[k] = tt / *HH(k, k);
421:   }

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

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

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

449:   KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP);

451:   /* add solution to previous solution */
452:   /* put updated solution into vdest.*/
453:   VecCopy(vguess, vdest);
454:   VecAXPY(vdest, 1.0, VEC_TEMP);
455:   return 0;
456: }

458: /*
459:     KSPLGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
460:                             Return new residual.

462:     input parameters:

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

469:     output parameters:
470: .        res - the new residual

472:  */
473: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
474: {
475:   PetscScalar *hh, *cc, *ss, tt;
476:   PetscInt     j;
477:   KSP_LGMRES  *lgmres = (KSP_LGMRES *)(ksp->data);

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

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

488:   for (j = 1; j <= it; j++) {
489:     tt  = *hh;
490:     *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
491:     hh++;
492:     *hh = *cc++ * *hh - (*ss++ * tt);
493:     /* hh, cc, and ss have all been incremented one by end of loop */
494:   }

496:   /*
497:     compute the new plane rotation, and apply it to:
498:      1) the right-hand-side of the Hessenberg system (GRS)
499:         note: it affects GRS(it) and GRS(it+1)
500:      2) the new column of the Hessenberg matrix
501:         note: it affects HH(it,it) which is currently pointed to
502:         by hh and HH(it+1, it) (*(hh+1))
503:     thus obtaining the updated value of the residual...
504:   */

506:   /* compute new plane rotation */

508:   if (!hapend) {
509:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
510:     if (tt == 0.0) {
511:       ksp->reason = KSP_DIVERGED_NULL;
512:       return 0;
513:     }
514:     *cc = *hh / tt;       /* new cosine value */
515:     *ss = *(hh + 1) / tt; /* new sine value */

517:     /* apply to 1) and 2) */
518:     *GRS(it + 1) = -(*ss * *GRS(it));
519:     *GRS(it)     = PetscConj(*cc) * *GRS(it);
520:     *hh          = PetscConj(*cc) * *hh + *ss * *(hh + 1);

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

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

532:     *res = 0.0;
533:   }
534:   return 0;
535: }

537: /*

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

542: */
543: static PetscErrorCode KSPLGMRESGetNewVectors(KSP ksp, PetscInt it)
544: {
545:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
546:   PetscInt    nwork  = lgmres->nwork_alloc; /* number of work vector chunks allocated */
547:   PetscInt    nalloc;                       /* number to allocate */
548:   PetscInt    k;

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

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

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

560:   /* work vectors */
561:   KSPCreateVecs(ksp, nalloc, &lgmres->user_work[nwork], 0, NULL);
562:   /* specify size of chunk allocated */
563:   lgmres->mwork_alloc[nwork] = nalloc;

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

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

569:   /* increment the number of work vector chunks */
570:   lgmres->nwork_alloc++;
571:   return 0;
572: }

574: /*
575:    KSPBuildSolution_LGMRES

577:      Input Parameter:
578: .     ksp - the Krylov space object
579: .     ptr-

581:    Output Parameter:
582: .     result - the solution

584:    Note: this calls KSPLGMRESBuildSoln - the same function that KSPLGMRESCycle
585:    calls directly.

587: */
588: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp, Vec ptr, Vec *result)
589: {
590:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

592:   if (!ptr) {
593:     if (!lgmres->sol_temp) { VecDuplicate(ksp->vec_sol, &lgmres->sol_temp); }
594:     ptr = lgmres->sol_temp;
595:   }
596:   if (!lgmres->nrs) {
597:     /* allocate the work area */
598:     PetscMalloc1(lgmres->max_k, &lgmres->nrs);
599:   }

601:   KSPLGMRESBuildSoln(lgmres->nrs, ksp->vec_sol, ptr, ksp, lgmres->it);
602:   if (result) *result = ptr;
603:   return 0;
604: }

606: PetscErrorCode KSPView_LGMRES(KSP ksp, PetscViewer viewer)
607: {
608:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
609:   PetscBool   iascii;

611:   KSPView_GMRES(ksp, viewer);
612:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
613:   if (iascii) {
614:     /*LGMRES_MOD */
615:     PetscViewerASCIIPrintf(viewer, "  aug. dimension=%" PetscInt_FMT "\n", lgmres->aug_dim);
616:     if (lgmres->approx_constant) PetscViewerASCIIPrintf(viewer, "  approx. space size was kept constant.\n");
617:     PetscViewerASCIIPrintf(viewer, "  number of matvecs=%" PetscInt_FMT "\n", lgmres->matvecs);
618:   }
619:   return 0;
620: }

622: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
623: {
624:   PetscInt    aug;
625:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
626:   PetscBool   flg    = PETSC_FALSE;

628:   KSPSetFromOptions_GMRES(ksp, PetscOptionsObject);
629:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP LGMRES Options");
630:   PetscOptionsBool("-ksp_lgmres_constant", "Use constant approx. space size", "KSPGMRESSetConstant", lgmres->approx_constant, &lgmres->approx_constant, NULL);
631:   PetscOptionsInt("-ksp_lgmres_augment", "Number of error approximations to augment the Krylov space with", "KSPLGMRESSetAugDim", lgmres->aug_dim, &aug, &flg);
632:   if (flg) KSPLGMRESSetAugDim(ksp, aug);
633:   PetscOptionsHeadEnd();
634:   return 0;
635: }

637: /*functions for extra lgmres options here*/
638: static PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)
639: {
640:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

642:   lgmres->approx_constant = PETSC_TRUE;
643:   return 0;
644: }

646: static PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp, PetscInt aug_dim)
647: {
648:   KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;

652:   lgmres->aug_dim = aug_dim;
653:   return 0;
654: }

656: /*MC
657:     KSPLGMRES - Augments the standard GMRES approximation space with approximations to
658:                 the error from previous restart cycles.

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

673:     To run LGMRES(m, k) as described in [1], use:
674:        -ksp_gmres_restart <m+k>
675:        -ksp_lgmres_augment <k>

677:   Level: beginner

679:    Note:
680:     Supports both left and right preconditioning, but not symmetric.

682:   Developer Note:
683:     This object is subclassed off of `KSPGMRES`

685:   Contributed by:
686:   Allison Baker

688:    References:
689: .  [1] - 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).

691: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPGMRES`,
692:           `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
693:           `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
694:           `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPLGMRESSetAugDim()`,
695:           `KSPGMRESSetConstant()`
696: M*/

698: PETSC_EXTERN PetscErrorCode KSPCreate_LGMRES(KSP ksp)
699: {
700:   KSP_LGMRES *lgmres;

702:   PetscNew(&lgmres);

704:   ksp->data               = (void *)lgmres;
705:   ksp->ops->buildsolution = KSPBuildSolution_LGMRES;

707:   ksp->ops->setup                        = KSPSetUp_LGMRES;
708:   ksp->ops->solve                        = KSPSolve_LGMRES;
709:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
710:   ksp->ops->view                         = KSPView_LGMRES;
711:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
712:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
713:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

715:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
716:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
717:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

719:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES);
720:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES);
721:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetOrthogonalization_C", KSPGMRESGetOrthogonalization_GMRES);
722:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES);
723:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_GMRES);
724:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES);
725:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES);
726:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetCGSRefinementType_C", KSPGMRESGetCGSRefinementType_GMRES);

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

732:   /*defaults */
733:   lgmres->haptol         = 1.0e-30;
734:   lgmres->q_preallocate  = 0;
735:   lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
736:   lgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
737:   lgmres->nrs            = NULL;
738:   lgmres->sol_temp       = NULL;
739:   lgmres->max_k          = LGMRES_DEFAULT_MAXK;
740:   lgmres->Rsvd           = NULL;
741:   lgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
742:   lgmres->orthogwork     = NULL;

744:   /*LGMRES_MOD - new defaults */
745:   lgmres->aug_dim         = LGMRES_DEFAULT_AUGDIM;
746:   lgmres->aug_ct          = 0; /* start with no aug vectors */
747:   lgmres->approx_constant = PETSC_FALSE;
748:   lgmres->matvecs         = 0;
749:   return 0;
750: }