Actual source code: pipefgmres.c

  1: /*
  2:   Contributed by Patrick Sanan and Sascha M. Schnepp
  3: */

  5: #include <../src/ksp/ksp/impls/gmres/pipefgmres/pipefgmresimpl.h>

  7: static PetscBool  cited      = PETSC_FALSE;
  8: static const char citation[] = "@article{SSM2016,\n"
  9:                                "  author = {P. Sanan and S.M. Schnepp and D.A. May},\n"
 10:                                "  title = {Pipelined, Flexible Krylov Subspace Methods},\n"
 11:                                "  journal = {SIAM Journal on Scientific Computing},\n"
 12:                                "  volume = {38},\n"
 13:                                "  number = {5},\n"
 14:                                "  pages = {C441-C470},\n"
 15:                                "  year = {2016},\n"
 16:                                "  doi = {10.1137/15M1049130},\n"
 17:                                "  URL = {http://dx.doi.org/10.1137/15M1049130},\n"
 18:                                "  eprint = {http://dx.doi.org/10.1137/15M1049130}\n"
 19:                                "}\n";

 21: #define PIPEFGMRES_DELTA_DIRECTIONS 10
 22: #define PIPEFGMRES_DEFAULT_MAXK     30

 24: static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP, PetscInt);
 25: static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP, PetscInt, PetscBool *, PetscReal *);
 26: static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);
 27: extern PetscErrorCode KSPReset_PIPEFGMRES(KSP);

 29: /*

 31:     KSPSetUp_PIPEFGMRES - Sets up the workspace needed by pipefgmres.

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

 36: */
 37: static PetscErrorCode KSPSetUp_PIPEFGMRES(KSP ksp)
 38: {
 39:   PetscInt        k;
 40:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
 41:   const PetscInt  max_k      = pipefgmres->max_k;

 43:   KSPSetUp_GMRES(ksp);

 45:   PetscMalloc1((VEC_OFFSET + max_k), &pipefgmres->prevecs);
 46:   PetscMalloc1((VEC_OFFSET + max_k), &pipefgmres->prevecs_user_work);

 48:   KSPCreateVecs(ksp, pipefgmres->vv_allocated, &pipefgmres->prevecs_user_work[0], 0, NULL);
 49:   for (k = 0; k < pipefgmres->vv_allocated; k++) pipefgmres->prevecs[k] = pipefgmres->prevecs_user_work[0][k];

 51:   PetscMalloc1((VEC_OFFSET + max_k), &pipefgmres->zvecs);
 52:   PetscMalloc1((VEC_OFFSET + max_k), &pipefgmres->zvecs_user_work);

 54:   PetscMalloc1((VEC_OFFSET + max_k), &pipefgmres->redux);

 56:   KSPCreateVecs(ksp, pipefgmres->vv_allocated, &pipefgmres->zvecs_user_work[0], 0, NULL);
 57:   for (k = 0; k < pipefgmres->vv_allocated; k++) pipefgmres->zvecs[k] = pipefgmres->zvecs_user_work[0][k];

 59:   return 0;
 60: }

 62: /*

 64:     KSPPIPEFGMRESCycle - Run pipefgmres, possibly with restart.  Return residual
 65:                   history if requested.

 67:     input parameters:
 68: .        pipefgmres  - structure containing parameters and work areas

 70:     output parameters:
 71: .        itcount - number of iterations used.  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: static PetscErrorCode KSPPIPEFGMRESCycle(PetscInt *itcount, KSP ksp)
 80: {
 81:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)(ksp->data);
 82:   PetscReal       res_norm;
 83:   PetscReal       hapbnd, tt;
 84:   PetscScalar    *hh, *hes, *lhh, shift = pipefgmres->shift;
 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 = pipefgmres->max_k; /* max # of directions Krylov space */
 88:   PetscInt        i, j, k;
 89:   Mat             Amat, Pmat;
 90:   Vec             Q, W;                      /* Pipelining vectors */
 91:   Vec            *redux = pipefgmres->redux; /* workspace for single reduction */

 93:   if (itcount) *itcount = 0;

 95:   /* Assign simpler names to these vectors, allocated as pipelining workspace */
 96:   Q = VEC_Q;
 97:   W = VEC_W;

 99:   /* Allocate memory for orthogonalization work (freed in the GMRES Destroy routine)*/
100:   /* Note that we add an extra value here to allow for a single reduction */
101:   if (!pipefgmres->orthogwork) PetscMalloc1(pipefgmres->max_k + 2, &pipefgmres->orthogwork);
102:   lhh = pipefgmres->orthogwork;

104:   /* Number of pseudo iterations since last restart is the number
105:      of prestart directions */
106:   loc_it = 0;

108:   /* note: (pipefgmres->it) is always set one less than (loc_it) It is used in
109:      KSPBUILDSolution_PIPEFGMRES, where it is passed to KSPPIPEFGMRESBuildSoln.
110:      Note that when KSPPIPEFGMRESBuildSoln is called from this function,
111:      (loc_it -1) is passed, so the two are equivalent */
112:   pipefgmres->it = (loc_it - 1);

114:   /* initial residual is in VEC_VV(0)  - compute its norm*/
115:   VecNorm(VEC_VV(0), NORM_2, &res_norm);

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

121:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
122:   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res_norm;
123:   else ksp->rnorm = 0;
124:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
125:   KSPLogResidualHistory(ksp, ksp->rnorm);
126:   KSPMonitor(ksp, ksp->its, ksp->rnorm);

128:   /* check for the convergence - maybe the current guess is good enough */
129:   (*ksp->converged)(ksp, ksp->its, ksp->rnorm, &ksp->reason, ksp->cnvP);
130:   if (ksp->reason) {
131:     if (itcount) *itcount = 0;
132:     return 0;
133:   }

135:   /* scale VEC_VV (the initial residual) */
136:   VecScale(VEC_VV(0), 1.0 / res_norm);

138:   /* Fill the pipeline */
139:   KSP_PCApply(ksp, VEC_VV(loc_it), PREVEC(loc_it));
140:   PCGetOperators(ksp->pc, &Amat, &Pmat);
141:   KSP_MatMult(ksp, Amat, PREVEC(loc_it), ZVEC(loc_it));
142:   VecAXPY(ZVEC(loc_it), -shift, VEC_VV(loc_it)); /* Note shift */

144:   /* 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:   while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
148:     if (loc_it) {
149:       KSPLogResidualHistory(ksp, res_norm);
150:       KSPMonitor(ksp, ksp->its, res_norm);
151:     }
152:     pipefgmres->it = (loc_it - 1);

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

161:     /* Note that these inner products are with "Z" now, so
162:        in particular, lhh[loc_it] is the 'barred' or 'shifted' value,
163:        not the value from the equivalent FGMRES run (even in exact arithmetic)
164:        That is, the H we need for the Arnoldi relation is different from the
165:        coefficients we use in the orthogonalization process,because of the shift */

167:     /* Do some local twiddling to allow for a single reduction */
168:     for (i = 0; i < loc_it + 1; i++) redux[i] = VEC_VV(i);
169:     redux[loc_it + 1] = ZVEC(loc_it);

171:     /* note the extra dot product which ends up in lh[loc_it+1], which computes ||z||^2 */
172:     VecMDotBegin(ZVEC(loc_it), loc_it + 2, redux, lhh);

174:     /* Start the split reduction (This actually calls the MPI_Iallreduce, otherwise, the reduction is simply delayed until the "end" call)*/
175:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)ZVEC(loc_it)));

177:     /* The work to be overlapped with the inner products follows.
178:        This is application of the preconditioner and the operator
179:        to compute intermediate quantites which will be combined (locally)
180:        with the results of the inner products.
181:        */
182:     KSP_PCApply(ksp, ZVEC(loc_it), Q);
183:     PCGetOperators(ksp->pc, &Amat, &Pmat);
184:     KSP_MatMult(ksp, Amat, Q, W);

186:     /* Compute inner products of the new direction with previous directions,
187:        and the norm of the to-be-orthogonalized direction "Z".
188:        This information is enough to build the required entries
189:        of H. The inner product with VEC_VV(it_loc) is
190:        *different* than in the standard FGMRES and need to be dealt with specially.
191:        That is, for standard FGMRES the orthogonalization coefficients are the same
192:        as the coefficients used in the Arnoldi relation to reconstruct, but here this
193:        is not true (albeit only for the one entry of H which we "unshift" below. */

195:     /* Finish the dot product, retrieving the extra entry */
196:     VecMDotEnd(ZVEC(loc_it), loc_it + 2, redux, lhh);
197:     tt = PetscRealPart(lhh[loc_it + 1]);

199:     /* Hessenberg entries, and entries for (naive) classical Graham-Schmidt
200:       Note that the Hessenberg entries require a shift, as these are for the
201:       relation AU = VH, which is wrt unshifted basis vectors */
202:     hh  = HH(0, loc_it);
203:     hes = HES(0, loc_it);
204:     for (j = 0; j < loc_it; j++) {
205:       hh[j]  = lhh[j];
206:       hes[j] = lhh[j];
207:     }
208:     hh[loc_it]  = lhh[loc_it] + shift;
209:     hes[loc_it] = lhh[loc_it] + shift;

211:     /* we delay applying the shift here */
212:     for (j = 0; j <= loc_it; j++) { lhh[j] = -lhh[j]; /* flip sign */ }

214:     /* Compute the norm of the un-normalized new direction using the rearranged formula
215:        Note that these are shifted ("barred") quantities */
216:     for (k = 0; k <= loc_it; k++) tt -= ((PetscReal)(PetscAbsScalar(lhh[k]) * PetscAbsScalar(lhh[k])));
217:     /* On AVX512 this is accumulating roundoff errors for eg: tt=-2.22045e-16 */
218:     if ((tt < 0.0) && tt > -PETSC_SMALL) tt = 0.0;
219:     if (tt < 0.0) {
220:       /* If we detect square root breakdown in the norm, we must restart the algorithm.
221:          Here this means we simply break the current loop and reconstruct the solution
222:          using the basis we have computed thus far. Note that by breaking immediately,
223:          we do not update the iteration count, so computation done in this iteration
224:          should be disregarded.
225:          */
226:       PetscInfo(ksp, "Restart due to square root breakdown at it = %" PetscInt_FMT ", tt=%g\n", ksp->its, (double)tt);
227:       break;
228:     } else {
229:       tt = PetscSqrtReal(tt);
230:     }

232:     /* new entry in hessenburg is the 2-norm of our new direction */
233:     hh[loc_it + 1]  = tt;
234:     hes[loc_it + 1] = tt;

236:     /* The recurred computation for the new direction
237:        The division by tt is delayed to the happy breakdown check later
238:        Note placement BEFORE the unshift
239:        */
240:     VecCopy(ZVEC(loc_it), VEC_VV(loc_it + 1));
241:     VecMAXPY(VEC_VV(loc_it + 1), loc_it + 1, lhh, &VEC_VV(0));
242:     /* (VEC_VV(loc_it+1) is not normalized yet) */

244:     /* The recurred computation for the preconditioned vector (u) */
245:     VecCopy(Q, PREVEC(loc_it + 1));
246:     VecMAXPY(PREVEC(loc_it + 1), loc_it + 1, lhh, &PREVEC(0));
247:     VecScale(PREVEC(loc_it + 1), 1.0 / tt);

249:     /* Unshift an entry in the GS coefficients ("removing the bar") */
250:     lhh[loc_it] -= shift;

252:     /* The recurred computation for z (Au)
253:        Note placement AFTER the "unshift" */
254:     VecCopy(W, ZVEC(loc_it + 1));
255:     VecMAXPY(ZVEC(loc_it + 1), loc_it + 1, lhh, &ZVEC(0));
256:     VecScale(ZVEC(loc_it + 1), 1.0 / tt);

258:     /* Happy Breakdown Check */
259:     hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
260:     /* RS(loc_it) contains the res_norm from the last iteration  */
261:     hapbnd = PetscMin(pipefgmres->haptol, hapbnd);
262:     if (tt > hapbnd) {
263:       /* scale new direction by its norm  */
264:       VecScale(VEC_VV(loc_it + 1), 1.0 / tt);
265:     } else {
266:       /* This happens when the solution is exactly reached. */
267:       /* So there is no new direction... */
268:       VecSet(VEC_TEMP, 0.0); /* set VEC_TEMP to 0 */
269:       hapend = PETSC_TRUE;
270:     }
271:     /* note that for pipefgmres we could get HES(loc_it+1, loc_it)  = 0 and the
272:        current solution would not be exact if HES was singular.  Note that
273:        HH non-singular implies that HES is not singular, and HES is guaranteed
274:        to be nonsingular when PREVECS are linearly independent and A is
275:        nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
276:        of HES). So we should really add a check to verify that HES is nonsingular.*/

278:     /* Note that to be thorough, in debug mode, one could call a LAPACK routine
279:        here to check that the hessenberg matrix is indeed non-singular (since
280:        FGMRES does not guarantee this) */

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

287:     loc_it++;
288:     pipefgmres->it = (loc_it - 1); /* Add this here in case it has converged */

290:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
291:     ksp->its++;
292:     if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res_norm;
293:     else ksp->rnorm = 0;
294:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

296:     (*ksp->converged)(ksp, ksp->its, ksp->rnorm, &ksp->reason, ksp->cnvP);

298:     /* Catch error in happy breakdown and signal convergence and break from loop */
299:     if (hapend) {
300:       if (!ksp->reason) {
302:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
303:         break;
304:       }
305:     }
306:   }
307:   /* END OF ITERATION LOOP */
308:   KSPLogResidualHistory(ksp, ksp->rnorm);

310:   /*
311:      Monitor if we know that we will not return for a restart */
312:   if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) KSPMonitor(ksp, ksp->its, ksp->rnorm);

314:   if (itcount) *itcount = loc_it;

316:   /*
317:     Down here we have to solve for the "best" coefficients of the Krylov
318:     columns, add the solution values together, and possibly unwind the
319:     preconditioning from the solution
320:    */

322:   /* Form the solution (or the solution so far) */
323:   /* Note: must pass in (loc_it-1) for iteration count so that KSPPIPEGMRESIIBuildSoln
324:      properly navigates */

326:   KSPPIPEFGMRESBuildSoln(RS(0), ksp->vec_sol, ksp->vec_sol, ksp, loc_it - 1);

328:   return 0;
329: }

331: /*
332:     KSPSolve_PIPEFGMRES - This routine applies the PIPEFGMRES method.

334:    Input Parameter:
335: .     ksp - the Krylov space object that was set to use pipefgmres

337:    Output Parameter:
338: .     outits - number of iterations used

340: */
341: static PetscErrorCode KSPSolve_PIPEFGMRES(KSP ksp)
342: {
343:   PetscInt        its, itcount;
344:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
345:   PetscBool       guess_zero = ksp->guess_zero;

347:   /* We have not checked these routines for use with complex numbers. The inner products
348:      are likely not defined correctly for that case */

351:   PetscCitationsRegister(citation, &cited);

354:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
355:   ksp->its = 0;
356:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

358:   itcount     = 0;
359:   ksp->reason = KSP_CONVERGED_ITERATING;
360:   while (!ksp->reason) {
361:     KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs);
362:     KSPPIPEFGMRESCycle(&its, ksp);
363:     itcount += its;
364:     if (itcount >= ksp->max_it) {
365:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
366:       break;
367:     }
368:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
369:   }
370:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
371:   return 0;
372: }

374: static PetscErrorCode KSPDestroy_PIPEFGMRES(KSP ksp)
375: {
376:   KSPReset_PIPEFGMRES(ksp);
377:   KSPDestroy_GMRES(ksp);
378:   return 0;
379: }

381: /*
382:     KSPPIPEFGMRESBuildSoln - create the solution from the starting vector and the
383:                       current iterates.

385:     Input parameters:
386:         nrs - work area of size it + 1.
387:         vguess  - index of initial guess
388:         vdest - index of result.  Note that vguess may == vdest (replace
389:                 guess with the solution).
390:         it - HH upper triangular part is a block of size (it+1) x (it+1)

392:      This is an internal routine that knows about the PIPEFGMRES internals.
393:  */
394: static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar *nrs, Vec vguess, Vec vdest, KSP ksp, PetscInt it)
395: {
396:   PetscScalar     tt;
397:   PetscInt        k, j;
398:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)(ksp->data);

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

402:   if (it < 0) {                        /* no pipefgmres steps have been performed */
403:     VecCopy(vguess, vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
404:     return 0;
405:   }

407:   /* solve the upper triangular system - RS is the right side and HH is
408:      the upper triangular matrix  - put soln in nrs */
409:   if (*HH(it, it) != 0.0) nrs[it] = *RS(it) / *HH(it, it);
410:   else nrs[it] = 0.0;

412:   for (k = it - 1; k >= 0; k--) {
413:     tt = *RS(k);
414:     for (j = k + 1; j <= it; j++) tt -= *HH(k, j) * nrs[j];
415:     nrs[k] = tt / *HH(k, k);
416:   }

418:   /* Accumulate the correction to the solution of the preconditioned problem in VEC_TEMP */
419:   VecZeroEntries(VEC_TEMP);
420:   VecMAXPY(VEC_TEMP, it + 1, nrs, &PREVEC(0));

422:   /* add solution to previous solution */
423:   if (vdest == vguess) {
424:     VecAXPY(vdest, 1.0, VEC_TEMP);
425:   } else {
426:     VecWAXPY(vdest, 1.0, VEC_TEMP, vguess);
427:   }
428:   return 0;
429: }

431: /*

433:     KSPPIPEFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
434:                             Return new residual.

436:     input parameters:

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

443:     output parameters:
444: .        res - the new residual

446:  */
447: /*
448: .  it - column of the Hessenberg that is complete, PIPEFGMRES is actually computing two columns ahead of this
449:  */
450: static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool *hapend, PetscReal *res)
451: {
452:   PetscScalar    *hh, *cc, *ss, *rs;
453:   PetscInt        j;
454:   PetscReal       hapbnd;
455:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)(ksp->data);

457:   hh = HH(0, it); /* pointer to beginning of column to update */
458:   cc = CC(0);     /* beginning of cosine rotations */
459:   ss = SS(0);     /* beginning of sine rotations */
460:   rs = RS(0);     /* right hand side of least squares system */

462:   /* The Hessenberg matrix is now correct through column it, save that form for possible spectral analysis */
463:   for (j = 0; j <= it + 1; j++) *HES(j, it) = hh[j];

465:   /* check for the happy breakdown */
466:   hapbnd = PetscMin(PetscAbsScalar(hh[it + 1] / rs[it]), pipefgmres->haptol);
467:   if (PetscAbsScalar(hh[it + 1]) < hapbnd) {
468:     PetscInfo(ksp, "Detected happy breakdown, current hapbnd = %14.12e H(%" PetscInt_FMT ",%" PetscInt_FMT ") = %14.12e\n", (double)hapbnd, it + 1, it, (double)PetscAbsScalar(*HH(it + 1, it)));
469:     *hapend = PETSC_TRUE;
470:   }

472:   /* Apply all the previously computed plane rotations to the new column
473:      of the Hessenberg matrix */
474:   /* Note: this uses the rotation [conj(c)  s ; -s   c], c= cos(theta), s= sin(theta),
475:      and some refs have [c   s ; -conj(s)  c] (don't be confused!) */

477:   for (j = 0; j < it; j++) {
478:     PetscScalar hhj = hh[j];
479:     hh[j]           = PetscConj(cc[j]) * hhj + ss[j] * hh[j + 1];
480:     hh[j + 1]       = -ss[j] * hhj + cc[j] * hh[j + 1];
481:   }

483:   /*
484:     compute the new plane rotation, and apply it to:
485:      1) the right-hand-side of the Hessenberg system (RS)
486:         note: it affects RS(it) and RS(it+1)
487:      2) the new column of the Hessenberg matrix
488:         note: it affects HH(it,it) which is currently pointed to
489:         by hh and HH(it+1, it) (*(hh+1))
490:     thus obtaining the updated value of the residual...
491:   */

493:   /* compute new plane rotation */

495:   if (!*hapend) {
496:     PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it + 1])));
497:     if (delta == 0.0) {
498:       ksp->reason = KSP_DIVERGED_NULL;
499:       return 0;
500:     }

502:     cc[it] = hh[it] / delta;     /* new cosine value */
503:     ss[it] = hh[it + 1] / delta; /* new sine value */

505:     hh[it]     = PetscConj(cc[it]) * hh[it] + ss[it] * hh[it + 1];
506:     rs[it + 1] = -ss[it] * rs[it];
507:     rs[it]     = PetscConj(cc[it]) * rs[it];
508:     *res       = PetscAbsScalar(rs[it + 1]);
509:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
510:             another rotation matrix (so RH doesn't change).  The new residual is
511:             always the new sine term times the residual from last time (RS(it)),
512:             but now the new sine rotation would be zero...so the residual should
513:             be zero...so we will multiply "zero" by the last residual.  This might
514:             not be exactly what we want to do here -could just return "zero". */

516:     *res = 0.0;
517:   }
518:   return 0;
519: }

521: /*
522:    KSPBuildSolution_PIPEFGMRES

524:      Input Parameter:
525: .     ksp - the Krylov space object
526: .     ptr-

528:    Output Parameter:
529: .     result - the solution

531:    Note: this calls KSPPIPEFGMRESBuildSoln - the same function that KSPPIPEFGMRESCycle
532:    calls directly.

534: */
535: PetscErrorCode KSPBuildSolution_PIPEFGMRES(KSP ksp, Vec ptr, Vec *result)
536: {
537:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;

539:   if (!ptr) {
540:     if (!pipefgmres->sol_temp) { VecDuplicate(ksp->vec_sol, &pipefgmres->sol_temp); }
541:     ptr = pipefgmres->sol_temp;
542:   }
543:   if (!pipefgmres->nrs) {
544:     /* allocate the work area */
545:     PetscMalloc1(pipefgmres->max_k, &pipefgmres->nrs);
546:   }

548:   KSPPIPEFGMRESBuildSoln(pipefgmres->nrs, ksp->vec_sol, ptr, ksp, pipefgmres->it);
549:   if (result) *result = ptr;
550:   return 0;
551: }

553: PetscErrorCode KSPSetFromOptions_PIPEFGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
554: {
555:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
556:   PetscBool       flg;
557:   PetscScalar     shift;

559:   KSPSetFromOptions_GMRES(ksp, PetscOptionsObject);
560:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP pipelined FGMRES Options");
561:   PetscOptionsScalar("-ksp_pipefgmres_shift", "shift parameter", "KSPPIPEFGMRESSetShift", pipefgmres->shift, &shift, &flg);
562:   if (flg) KSPPIPEFGMRESSetShift(ksp, shift);
563:   PetscOptionsHeadEnd();
564:   return 0;
565: }

567: PetscErrorCode KSPView_PIPEFGMRES(KSP ksp, PetscViewer viewer)
568: {
569:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
570:   PetscBool       iascii, isstring;

572:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
573:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring);

575:   if (iascii) {
576:     PetscViewerASCIIPrintf(viewer, "  restart=%" PetscInt_FMT "\n", pipefgmres->max_k);
577:     PetscViewerASCIIPrintf(viewer, "  happy breakdown tolerance %g\n", (double)pipefgmres->haptol);
578: #if defined(PETSC_USE_COMPLEX)
579:     PetscViewerASCIIPrintf(viewer, "  shift=%g+%gi\n", (double)PetscRealPart(pipefgmres->shift), (double)PetscImaginaryPart(pipefgmres->shift));
580: #else
581:     PetscViewerASCIIPrintf(viewer, "  shift=%g\n", (double)pipefgmres->shift);
582: #endif
583:   } else if (isstring) {
584:     PetscViewerStringSPrintf(viewer, "restart %" PetscInt_FMT, pipefgmres->max_k);
585: #if defined(PETSC_USE_COMPLEX)
586:     PetscViewerStringSPrintf(viewer, "   shift=%g+%gi\n", (double)PetscRealPart(pipefgmres->shift), (double)PetscImaginaryPart(pipefgmres->shift));
587: #else
588:     PetscViewerStringSPrintf(viewer, "   shift=%g\n", (double)pipefgmres->shift);
589: #endif
590:   }
591:   return 0;
592: }

594: PetscErrorCode KSPReset_PIPEFGMRES(KSP ksp)
595: {
596:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
597:   PetscInt        i;

599:   PetscFree(pipefgmres->prevecs);
600:   PetscFree(pipefgmres->zvecs);
601:   for (i = 0; i < pipefgmres->nwork_alloc; i++) {
602:     VecDestroyVecs(pipefgmres->mwork_alloc[i], &pipefgmres->prevecs_user_work[i]);
603:     VecDestroyVecs(pipefgmres->mwork_alloc[i], &pipefgmres->zvecs_user_work[i]);
604:   }
605:   PetscFree(pipefgmres->prevecs_user_work);
606:   PetscFree(pipefgmres->zvecs_user_work);
607:   PetscFree(pipefgmres->redux);
608:   KSPReset_GMRES(ksp);
609:   return 0;
610: }

612: /*MC
613:    KSPPIPEFGMRES - Implements the Pipelined Generalized Minimal Residual method.

615:    A flexible, 1-stage pipelined variant of GMRES.

617:    Options Database Keys:
618: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
619: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
620: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
621: .   -ksp_pipefgmres_shift - the shift to use (defaults to 1. See KSPPIPEFGMRESSetShift()
622:                              vectors are allocated as needed)
623: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

625:    Level: intermediate

627:    Notes:

629:    This variant is not "explicitly normalized" like KSPPGMRES, and requires a shift parameter.

631:    A heuristic for choosing the shift parameter is the largest eigenvalue of the preconditioned operator.

633:    Only right preconditioning is supported (but this preconditioner may be nonlinear/variable/inexact, as with KSPFGMRES).
634:    MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
635:    See the FAQ on the PETSc website for details.

637:    Developer Notes:
638:     This class is subclassed off of KSPGMRES.

640:    Reference:
641:     P. Sanan, S.M. Schnepp, and D.A. May,
642:     "Pipelined, Flexible Krylov Subspace Methods,"
643:     SIAM Journal on Scientific Computing 2016 38:5, C441-C470,
644:     DOI: 10.1137/15M1049130

646: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPLGMRES`, `KSPPIPECG`, `KSPPIPECR`, `KSPPGMRES`, `KSPFGMRES`
647:           `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESMonitorKrylov()`, `KSPPIPEFGMRESSetShift()`
648: M*/

650: PETSC_EXTERN PetscErrorCode KSPCreate_PIPEFGMRES(KSP ksp)
651: {
652:   KSP_PIPEFGMRES *pipefgmres;

654:   PetscNew(&pipefgmres);

656:   ksp->data                              = (void *)pipefgmres;
657:   ksp->ops->buildsolution                = KSPBuildSolution_PIPEFGMRES;
658:   ksp->ops->setup                        = KSPSetUp_PIPEFGMRES;
659:   ksp->ops->solve                        = KSPSolve_PIPEFGMRES;
660:   ksp->ops->reset                        = KSPReset_PIPEFGMRES;
661:   ksp->ops->destroy                      = KSPDestroy_PIPEFGMRES;
662:   ksp->ops->view                         = KSPView_PIPEFGMRES;
663:   ksp->ops->setfromoptions               = KSPSetFromOptions_PIPEFGMRES;
664:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
665:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

667:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 3);
668:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

670:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES);
671:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES);
672:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESGetRestart_C", KSPGMRESGetRestart_GMRES);

674:   pipefgmres->nextra_vecs    = 1;
675:   pipefgmres->haptol         = 1.0e-30;
676:   pipefgmres->q_preallocate  = 0;
677:   pipefgmres->delta_allocate = PIPEFGMRES_DELTA_DIRECTIONS;
678:   pipefgmres->orthog         = NULL;
679:   pipefgmres->nrs            = NULL;
680:   pipefgmres->sol_temp       = NULL;
681:   pipefgmres->max_k          = PIPEFGMRES_DEFAULT_MAXK;
682:   pipefgmres->Rsvd           = NULL;
683:   pipefgmres->orthogwork     = NULL;
684:   pipefgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
685:   pipefgmres->shift          = 1.0;
686:   return 0;
687: }

689: static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP ksp, PetscInt it)
690: {
691:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;
692:   PetscInt        nwork      = pipefgmres->nwork_alloc; /* number of work vector chunks allocated */
693:   PetscInt        nalloc;                               /* number to allocate */
694:   PetscInt        k;

696:   nalloc = pipefgmres->delta_allocate; /* number of vectors to allocate
697:                                       in a single chunk */

699:   /* Adjust the number to allocate to make sure that we don't exceed the
700:      number of available slots (pipefgmres->vecs_allocated)*/
701:   if (it + VEC_OFFSET + nalloc >= pipefgmres->vecs_allocated) nalloc = pipefgmres->vecs_allocated - it - VEC_OFFSET;
702:   if (!nalloc) return 0;

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

706:   /* work vectors */
707:   KSPCreateVecs(ksp, nalloc, &pipefgmres->user_work[nwork], 0, NULL);
708:   for (k = 0; k < nalloc; k++) pipefgmres->vecs[it + VEC_OFFSET + k] = pipefgmres->user_work[nwork][k];
709:   /* specify size of chunk allocated */
710:   pipefgmres->mwork_alloc[nwork] = nalloc;

712:   /* preconditioned vectors (note we don't use VEC_OFFSET) */
713:   KSPCreateVecs(ksp, nalloc, &pipefgmres->prevecs_user_work[nwork], 0, NULL);
714:   for (k = 0; k < nalloc; k++) pipefgmres->prevecs[it + k] = pipefgmres->prevecs_user_work[nwork][k];

716:   KSPCreateVecs(ksp, nalloc, &pipefgmres->zvecs_user_work[nwork], 0, NULL);
717:   for (k = 0; k < nalloc; k++) pipefgmres->zvecs[it + k] = pipefgmres->zvecs_user_work[nwork][k];

719:   /* increment the number of work vector chunks */
720:   pipefgmres->nwork_alloc++;
721:   return 0;
722: }

724: /*@
725:   KSPPIPEFGMRESSetShift - Set the shift parameter for the flexible, pipelined GMRES solver.

727:   A heuristic is to set this to be comparable to the largest eigenvalue of the preconditioned operator. This can be acheived with PETSc itself by using a few iterations of a Krylov method. See KSPComputeEigenvalues (and note the caveats there).

729: Logically Collective on ksp

731: Input Parameters:
732: +  ksp - the Krylov space context
733: -  shift - the shift

735: Level: intermediate

737: Options Database:
738: . -ksp_pipefgmres_shift <shift> - set the shift parameter

740: .seealso: `KSPComputeEigenvalues()`
741: @*/
742: PetscErrorCode KSPPIPEFGMRESSetShift(KSP ksp, PetscScalar shift)
743: {
744:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES *)ksp->data;

748:   pipefgmres->shift = shift;
749:   return 0;
750: }