Actual source code: pgmres.c


  2: /*
  3:     This file implements PGMRES (a Pipelined Generalized Minimal Residual method)
  4: */

  6: #include <../src/ksp/ksp/impls/gmres/pgmres/pgmresimpl.h>
  7: #define PGMRES_DELTA_DIRECTIONS 10
  8: #define PGMRES_DEFAULT_MAXK     30

 10: static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP,PetscInt,PetscBool*,PetscReal*);
 11: static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 13: /*

 15:     KSPSetUp_PGMRES - Sets up the workspace needed by pgmres.

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

 20: */
 21: static PetscErrorCode KSPSetUp_PGMRES(KSP ksp)
 22: {
 23:   KSPSetUp_GMRES(ksp);
 24:   return 0;
 25: }

 27: /*

 29:     KSPPGMRESCycle - Run pgmres, possibly with restart.  Return residual
 30:                   history if requested.

 32:     input parameters:
 33: .        pgmres  - structure containing parameters and work areas

 35:     output parameters:
 36: .        itcount - number of iterations used.  If null, ignored.
 37: .        converged - 0 if not converged

 39:     Notes:
 40:     On entry, the value in vector VEC_VV(0) should be
 41:     the initial residual.

 43:  */
 44: static PetscErrorCode KSPPGMRESCycle(PetscInt *itcount,KSP ksp)
 45: {
 46:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);
 47:   PetscReal      res_norm,res,newnorm;
 48:   PetscInt       it     = 0,j,k;
 49:   PetscBool      hapend = PETSC_FALSE;

 51:   if (itcount) *itcount = 0;
 52:   VecNormalize(VEC_VV(0),&res_norm);
 53:   KSPCheckNorm(ksp,res_norm);
 54:   res    = res_norm;
 55:   *RS(0) = res_norm;

 57:   /* check for the convergence */
 58:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 59:   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
 60:   else ksp->rnorm = 0;
 61:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 62:   pgmres->it = it-2;
 63:   KSPLogResidualHistory(ksp,ksp->rnorm);
 64:   KSPMonitor(ksp,ksp->its,ksp->rnorm);
 65:   if (!res) {
 66:     ksp->reason = KSP_CONVERGED_ATOL;
 67:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
 68:     return 0;
 69:   }

 71:   (*ksp->converged)(ksp,ksp->its,ksp->rnorm,&ksp->reason,ksp->cnvP);
 72:   for (; !ksp->reason; it++) {
 73:     Vec Zcur,Znext;
 74:     if (pgmres->vv_allocated <= it + VEC_OFFSET + 1) {
 75:       KSPGMRESGetNewVectors(ksp,it+1);
 76:     }
 77:     /* VEC_VV(it-1) is orthogonal, it will be normalized once the VecNorm arrives. */
 78:     Zcur  = VEC_VV(it);         /* Zcur is not yet orthogonal, but the VecMDot to orthogonalize it has been started. */
 79:     Znext = VEC_VV(it+1);       /* This iteration will compute Znext, update with a deferred correction once we know how
 80:                                  * Zcur relates to the previous vectors, and start the reduction to orthogonalize it. */

 82:     if (it < pgmres->max_k+1 && ksp->its+1 < PetscMax(2,ksp->max_it)) { /* We don't know whether what we have computed is enough, so apply the matrix. */
 83:       KSP_PCApplyBAorAB(ksp,Zcur,Znext,VEC_TEMP_MATOP);
 84:     }

 86:     if (it > 1) {               /* Complete the pending reduction */
 87:       VecNormEnd(VEC_VV(it-1),NORM_2,&newnorm);
 88:       *HH(it-1,it-2) = newnorm;
 89:     }
 90:     if (it > 0) {               /* Finish the reduction computing the latest column of H */
 91:       VecMDotEnd(Zcur,it,&(VEC_VV(0)),HH(0,it-1));
 92:     }

 94:     if (it > 1) {
 95:       /* normalize the base vector from two iterations ago, basis is complete up to here */
 96:       VecScale(VEC_VV(it-1),1./ *HH(it-1,it-2));

 98:       KSPPGMRESUpdateHessenberg(ksp,it-2,&hapend,&res);
 99:       pgmres->it = it-2;
100:       ksp->its++;
101:       if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
102:       else ksp->rnorm = 0;

104:       (*ksp->converged)(ksp,ksp->its,ksp->rnorm,&ksp->reason,ksp->cnvP);
105:       if (it < pgmres->max_k+1 || ksp->reason || ksp->its == ksp->max_it) {  /* Monitor if we are done or still iterating, but not before a restart. */
106:         KSPLogResidualHistory(ksp,ksp->rnorm);
107:         KSPMonitor(ksp,ksp->its,ksp->rnorm);
108:       }
109:       if (ksp->reason) break;
110:       /* Catch error in happy breakdown and signal convergence and break from loop */
111:       if (hapend) {
113:         else {
114:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
115:           break;
116:         }
117:       }

119:       if (!(it < pgmres->max_k+1 && ksp->its < ksp->max_it)) break;

121:       /* The it-2 column of H was not scaled when we computed Zcur, apply correction */
122:       VecScale(Zcur,1./ *HH(it-1,it-2));
123:       /* And Znext computed in this iteration was computed using the under-scaled Zcur */
124:       VecScale(Znext,1./ *HH(it-1,it-2));

126:       /* In the previous iteration, we projected an unnormalized Zcur against the Krylov basis, so we need to fix the column of H resulting from that projection. */
127:       for (k=0; k<it; k++) *HH(k,it-1) /= *HH(it-1,it-2);
128:       /* When Zcur was projected against the Krylov basis, VV(it-1) was still not normalized, so fix that too. This
129:        * column is complete except for HH(it,it-1) which we won't know until the next iteration. */
130:       *HH(it-1,it-1) /= *HH(it-1,it-2);
131:     }

133:     if (it > 0) {
134:       PetscScalar *work;
135:       if (!pgmres->orthogwork) PetscMalloc1(pgmres->max_k + 2,&pgmres->orthogwork);
136:       work = pgmres->orthogwork;
137:       /* Apply correction computed by the VecMDot in the last iteration to Znext. The original form is
138:        *
139:        *   Znext -= sum_{j=0}^{i-1} Z[j+1] * H[j,i-1]
140:        *
141:        * where
142:        *
143:        *   Z[j] = sum_{k=0}^j V[k] * H[k,j-1]
144:        *
145:        * substituting
146:        *
147:        *   Znext -= sum_{j=0}^{i-1} sum_{k=0}^{j+1} V[k] * H[k,j] * H[j,i-1]
148:        *
149:        * rearranging the iteration space from row-column to column-row
150:        *
151:        *   Znext -= sum_{k=0}^i sum_{j=k-1}^{i-1} V[k] * H[k,j] * H[j,i-1]
152:        *
153:        * Note that column it-1 of HH is correct. For all previous columns, we must look at HES because HH has already
154:        * been transformed to upper triangular form.
155:        */
156:       for (k=0; k<it+1; k++) {
157:         work[k] = 0;
158:         for (j=PetscMax(0,k-1); j<it-1; j++) work[k] -= *HES(k,j) * *HH(j,it-1);
159:       }
160:       VecMAXPY(Znext,it+1,work,&VEC_VV(0));
161:       VecAXPY(Znext,-*HH(it-1,it-1),Zcur);

163:       /* Orthogonalize Zcur against existing basis vectors. */
164:       for (k=0; k<it; k++) work[k] = -*HH(k,it-1);
165:       VecMAXPY(Zcur,it,work,&VEC_VV(0));
166:       /* Zcur is now orthogonal, and will be referred to as VEC_VV(it) again, though it is still not normalized. */
167:       /* Begin computing the norm of the new vector, will be normalized after the MatMult in the next iteration. */
168:       VecNormBegin(VEC_VV(it),NORM_2,&newnorm);
169:     }

171:     /* Compute column of H (to the diagonal, but not the subdiagonal) to be able to orthogonalize the newest vector. */
172:     VecMDotBegin(Znext,it+1,&VEC_VV(0),HH(0,it));

174:     /* Start an asynchronous split-mode reduction, the result of the MDot and Norm will be collected on the next iteration. */
175:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Znext));
176:   }

178:   if (itcount) *itcount = it-1; /* Number of iterations actually completed. */

180:   /*
181:     Down here we have to solve for the "best" coefficients of the Krylov
182:     columns, add the solution values together, and possibly unwind the
183:     preconditioning from the solution
184:    */
185:   /* Form the solution (or the solution so far) */
186:   KSPPGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-2);
187:   return 0;
188: }

190: /*
191:     KSPSolve_PGMRES - This routine applies the PGMRES method.

193:    Input Parameter:
194: .     ksp - the Krylov space object that was set to use pgmres

196:    Output Parameter:
197: .     outits - number of iterations used

199: */
200: static PetscErrorCode KSPSolve_PGMRES(KSP ksp)
201: {
202:   PetscInt       its,itcount;
203:   KSP_PGMRES     *pgmres    = (KSP_PGMRES*)ksp->data;
204:   PetscBool      guess_zero = ksp->guess_zero;

207:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
208:   ksp->its = 0;
209:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

211:   itcount     = 0;
212:   ksp->reason = KSP_CONVERGED_ITERATING;
213:   while (!ksp->reason) {
214:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
215:     KSPPGMRESCycle(&its,ksp);
216:     itcount += its;
217:     if (itcount >= ksp->max_it) {
218:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
219:       break;
220:     }
221:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
222:   }
223:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
224:   return 0;
225: }

227: static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)
228: {
229:   KSPDestroy_GMRES(ksp);
230:   return 0;
231: }

233: /*
234:     KSPPGMRESBuildSoln - create the solution from the starting vector and the
235:                       current iterates.

237:     Input parameters:
238:         nrs - work area of size it + 1.
239:         vguess  - index of initial guess
240:         vdest - index of result.  Note that vguess may == vdest (replace
241:                 guess with the solution).
242:         it - HH upper triangular part is a block of size (it+1) x (it+1)

244:      This is an internal routine that knows about the PGMRES internals.
245:  */
246: static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
247: {
248:   PetscScalar    tt;
249:   PetscInt       k,j;
250:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);

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

254:   if (it < 0) {                                 /* no pgmres steps have been performed */
255:     VecCopy(vguess,vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
256:     return 0;
257:   }

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

264:   for (k=it-1; k>=0; k--) {
265:     tt = *RS(k);
266:     for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
267:     nrs[k] = tt / *HH(k,k);
268:   }

270:   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
271:   VecZeroEntries(VEC_TEMP);
272:   VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
273:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
274:   /* add solution to previous solution */
275:   if (vdest == vguess) {
276:     VecAXPY(vdest,1.0,VEC_TEMP);
277:   } else {
278:     VecWAXPY(vdest,1.0,VEC_TEMP,vguess);
279:   }
280:   return 0;
281: }

283: /*

285:     KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
286:                             Return new residual.

288:     input parameters:

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

295:     output parameters:
296: .        res - the new residual

298:  */
299: /*
300: .  it - column of the Hessenberg that is complete, PGMRES is actually computing two columns ahead of this
301:  */
302: static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res)
303: {
304:   PetscScalar    *hh,*cc,*ss,*rs;
305:   PetscInt       j;
306:   PetscReal      hapbnd;
307:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);

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

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

317:   /* check for the happy breakdown */
318:   hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pgmres->haptol);
319:   if (PetscAbsScalar(hh[it+1]) < hapbnd) {
320:     PetscInfo(ksp,"Detected happy breakdown, current hapbnd = %14.12e H(%D,%D) = %14.12e\n",(double)hapbnd,it+1,it,(double)PetscAbsScalar(*HH(it+1,it)));
321:     *hapend = PETSC_TRUE;
322:   }

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

329:   for (j=0; j<it; j++) {
330:     PetscScalar hhj = hh[j];
331:     hh[j]   = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
332:     hh[j+1] =          -ss[j] *hhj + cc[j]*hh[j+1];
333:   }

335:   /*
336:     compute the new plane rotation, and apply it to:
337:      1) the right-hand-side of the Hessenberg system (RS)
338:         note: it affects RS(it) and RS(it+1)
339:      2) the new column of the Hessenberg matrix
340:         note: it affects HH(it,it) which is currently pointed to
341:         by hh and HH(it+1, it) (*(hh+1))
342:     thus obtaining the updated value of the residual...
343:   */

345:   /* compute new plane rotation */

347:   if (!*hapend) {
348:     PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1])));
349:     if (delta == 0.0) {
350:       ksp->reason = KSP_DIVERGED_NULL;
351:       return 0;
352:     }

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

357:     hh[it]   = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1];
358:     rs[it+1] = -ss[it]*rs[it];
359:     rs[it]   = PetscConj(cc[it])*rs[it];
360:     *res     = PetscAbsScalar(rs[it+1]);
361:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
362:             another rotation matrix (so RH doesn't change).  The new residual is
363:             always the new sine term times the residual from last time (RS(it)),
364:             but now the new sine rotation would be zero...so the residual should
365:             be zero...so we will multiply "zero" by the last residual.  This might
366:             not be exactly what we want to do here -could just return "zero". */

368:     *res = 0.0;
369:   }
370:   return 0;
371: }

373: /*
374:    KSPBuildSolution_PGMRES

376:      Input Parameter:
377: .     ksp - the Krylov space object
378: .     ptr-

380:    Output Parameter:
381: .     result - the solution

383:    Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
384:    calls directly.

386: */
387: PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)
388: {
389:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)ksp->data;

391:   if (!ptr) {
392:     if (!pgmres->sol_temp) {
393:       VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);
394:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)pgmres->sol_temp);
395:     }
396:     ptr = pgmres->sol_temp;
397:   }
398:   if (!pgmres->nrs) {
399:     /* allocate the work area */
400:     PetscMalloc1(pgmres->max_k,&pgmres->nrs);
401:     PetscLogObjectMemory((PetscObject)ksp,pgmres->max_k*sizeof(PetscScalar));
402:   }

404:   KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);
405:   if (result) *result = ptr;
406:   return 0;
407: }

409: PetscErrorCode KSPSetFromOptions_PGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
410: {
411:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
412:   PetscOptionsHead(PetscOptionsObject,"KSP pipelined GMRES Options");
413:   PetscOptionsTail();
414:   return 0;
415: }

417: PetscErrorCode KSPReset_PGMRES(KSP ksp)
418: {
419:   KSPReset_GMRES(ksp);
420:   return 0;
421: }

423: /*MC
424:      KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.

426:    Options Database Keys:
427: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
428: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
429: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
430:                              vectors are allocated as needed)
431: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
432: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
433: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
434:                                    stability of the classical Gram-Schmidt  orthogonalization.
435: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

437:    Level: beginner

439:    Notes:
440:    MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
441:    See the FAQ on the PETSc website for details.

443:    Reference:
444:    Ghysels, Ashby, Meerbergen, Vanroose, Hiding global communication latencies in the GMRES algorithm on massively parallel machines, 2012.

446:    Developer Notes:
447:     This object is subclassed off of KSPGMRES

449: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES, KSPPIPECG, KSPPIPECR,
450:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
451:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
452:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov()
453: M*/

455: PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp)
456: {
457:   KSP_PGMRES     *pgmres;

459:   PetscNewLog(ksp,&pgmres);

461:   ksp->data                              = (void*)pgmres;
462:   ksp->ops->buildsolution                = KSPBuildSolution_PGMRES;
463:   ksp->ops->setup                        = KSPSetUp_PGMRES;
464:   ksp->ops->solve                        = KSPSolve_PGMRES;
465:   ksp->ops->reset                        = KSPReset_PGMRES;
466:   ksp->ops->destroy                      = KSPDestroy_PGMRES;
467:   ksp->ops->view                         = KSPView_GMRES;
468:   ksp->ops->setfromoptions               = KSPSetFromOptions_PGMRES;
469:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
470:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

472:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
473:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
474:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);

476:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
477:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
478:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
479:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
480:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
481:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
482:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

484:   pgmres->nextra_vecs    = 1;
485:   pgmres->haptol         = 1.0e-30;
486:   pgmres->q_preallocate  = 0;
487:   pgmres->delta_allocate = PGMRES_DELTA_DIRECTIONS;
488:   pgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
489:   pgmres->nrs            = NULL;
490:   pgmres->sol_temp       = NULL;
491:   pgmres->max_k          = PGMRES_DEFAULT_MAXK;
492:   pgmres->Rsvd           = NULL;
493:   pgmres->orthogwork     = NULL;
494:   pgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
495:   return 0;
496: }