Actual source code: pgmres.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: /*
  3:     This file implements PGMRES (a Pipelined Generalized Minimal Residual method)
  4: */

  6: #include <../src/ksp/ksp/impls/gmres/pgmres/pgmresimpl.h>       /*I  "petscksp.h"  I*/
  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: */
 23: static PetscErrorCode KSPSetUp_PGMRES(KSP ksp)
 24: {

 28:   KSPSetUp_GMRES(ksp);
 29:   return(0);
 30: }

 32: /*

 34:     KSPPGMRESCycle - Run pgmres, possibly with restart.  Return residual
 35:                   history if requested.

 37:     input parameters:
 38: .        pgmres  - structure containing parameters and work areas

 40:     output parameters:
 41: .        itcount - number of iterations used.  If null, ignored.
 42: .        converged - 0 if not converged

 44:     Notes:
 45:     On entry, the value in vector VEC_VV(0) should be
 46:     the initial residual.


 49:  */
 52: static PetscErrorCode KSPPGMRESCycle(PetscInt *itcount,KSP ksp)
 53: {
 54:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);
 55:   PetscReal      res_norm,res,newnorm;
 57:   PetscInt       it     = 0,j,k;
 58:   PetscBool      hapend = PETSC_FALSE;

 61:   if (itcount) *itcount = 0;
 62:   VecNormalize(VEC_VV(0),&res_norm);
 63:   KSPCheckNorm(ksp,res_norm);
 64:   res    = res_norm;
 65:   *RS(0) = res_norm;

 67:   /* check for the convergence */
 68:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 69:   ksp->rnorm = res;
 70:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 71:   pgmres->it = it-2;
 72:   KSPLogResidualHistory(ksp,res);
 73:   KSPMonitor(ksp,ksp->its,res);
 74:   if (!res) {
 75:     ksp->reason = KSP_CONVERGED_ATOL;
 76:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
 77:     return(0);
 78:   }

 80:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
 81:   for (; !ksp->reason; it++) {
 82:     Vec Zcur,Znext;
 83:     if (pgmres->vv_allocated <= it + VEC_OFFSET + 1) {
 84:       KSPGMRESGetNewVectors(ksp,it+1);
 85:     }
 86:     /* VEC_VV(it-1) is orthogonal, it will be normalized once the VecNorm arrives. */
 87:     Zcur  = VEC_VV(it);         /* Zcur is not yet orthogonal, but the VecMDot to orthogonalize it has been started. */
 88:     Znext = VEC_VV(it+1);       /* This iteration will compute Znext, update with a deferred correction once we know how
 89:                                  * Zcur relates to the previous vectors, and start the reduction to orthogonalize it. */

 91:     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. */
 92:       KSP_PCApplyBAorAB(ksp,Zcur,Znext,VEC_TEMP_MATOP);
 93:     }

 95:     if (it > 1) {               /* Complete the pending reduction */
 96:       VecNormEnd(VEC_VV(it-1),NORM_2,&newnorm);
 97:       *HH(it-1,it-2) = newnorm;
 98:     }
 99:     if (it > 0) {               /* Finish the reduction computing the latest column of H */
100:       VecMDotEnd(Zcur,it,&(VEC_VV(0)),HH(0,it-1));
101:     }

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

107:       KSPPGMRESUpdateHessenberg(ksp,it-2,&hapend,&res);
108:       pgmres->it = it-2;
109:       ksp->its++;
110:       ksp->rnorm = res;

112:       (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
113:       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. */
114:         KSPLogResidualHistory(ksp,res);
115:         KSPMonitor(ksp,ksp->its,res);
116:       }
117:       if (ksp->reason) break;
118:       /* Catch error in happy breakdown and signal convergence and break from loop */
119:       if (hapend) {
120:         if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",(double)res);
121:         else {
122:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
123:           break;
124:         }
125:       }

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

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

134:       /* 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. */
135:       for (k=0; k<it; k++) *HH(k,it-1) /= *HH(it-1,it-2);
136:       /* When Zcur was projected against the Krylov basis, VV(it-1) was still not normalized, so fix that too. This
137:        * column is complete except for HH(it,it-1) which we won't know until the next iteration. */
138:       *HH(it-1,it-1) /= *HH(it-1,it-2);
139:     }

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

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

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

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

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

188:   /*
189:     Down here we have to solve for the "best" coefficients of the Krylov
190:     columns, add the solution values together, and possibly unwind the
191:     preconditioning from the solution
192:    */
193:   /* Form the solution (or the solution so far) */
194:   KSPPGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-2);
195:   return(0);
196: }

198: /*
199:     KSPSolve_PGMRES - This routine applies the PGMRES method.


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

205:    Output Parameter:
206: .     outits - number of iterations used

208: */
211: static PetscErrorCode KSPSolve_PGMRES(KSP ksp)
212: {
214:   PetscInt       its,itcount;
215:   KSP_PGMRES     *pgmres    = (KSP_PGMRES*)ksp->data;
216:   PetscBool      guess_zero = ksp->guess_zero;

219:   if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
220:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
221:   ksp->its = 0;
222:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

224:   itcount     = 0;
225:   ksp->reason = KSP_CONVERGED_ITERATING;
226:   while (!ksp->reason) {
227:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
228:     KSPPGMRESCycle(&its,ksp);
229:     itcount += its;
230:     if (itcount >= ksp->max_it) {
231:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
232:       break;
233:     }
234:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
235:   }
236:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
237:   return(0);
238: }

242: static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)
243: {

247:   KSPDestroy_GMRES(ksp);
248:   return(0);
249: }

251: /*
252:     KSPPGMRESBuildSoln - create the solution from the starting vector and the
253:                       current iterates.

255:     Input parameters:
256:         nrs - work area of size it + 1.
257:         vguess  - index of initial guess
258:         vdest - index of result.  Note that vguess may == vdest (replace
259:                 guess with the solution).
260:         it - HH upper triangular part is a block of size (it+1) x (it+1)

262:      This is an internal routine that knows about the PGMRES internals.
263:  */
266: static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
267: {
268:   PetscScalar    tt;
270:   PetscInt       k,j;
271:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);

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

276:   if (it < 0) {                                 /* no pgmres steps have been performed */
277:     VecCopy(vguess,vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
278:     return(0);
279:   }

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

286:   for (k=it-1; k>=0; k--) {
287:     tt = *RS(k);
288:     for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
289:     nrs[k] = tt / *HH(k,k);
290:   }

292:   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
293:   VecZeroEntries(VEC_TEMP);
294:   VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
295:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
296:   /* add solution to previous solution */
297:   if (vdest == vguess) {
298:     VecAXPY(vdest,1.0,VEC_TEMP);
299:   } else {
300:     VecWAXPY(vdest,1.0,VEC_TEMP,vguess);
301:   }
302:   return(0);
303: }

305: /*

307:     KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
308:                             Return new residual.

310:     input parameters:

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

317:     output parameters:
318: .        res - the new residual

320:  */
323: /*
324: .  it - column of the Hessenberg that is complete, PGMRES is actually computing two columns ahead of this
325:  */
326: static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res)
327: {
328:   PetscScalar    *hh,*cc,*ss,*rs;
329:   PetscInt       j;
330:   PetscReal      hapbnd;
331:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);

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

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

343:   /* check for the happy breakdown */
344:   hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pgmres->haptol);
345:   if (PetscAbsScalar(hh[it+1]) < hapbnd) {
346:     PetscInfo4(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)));
347:     *hapend = PETSC_TRUE;
348:   }

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

355:   for (j=0; j<it; j++) {
356:     PetscScalar hhj = hh[j];
357:     hh[j]   = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
358:     hh[j+1] =          -ss[j] *hhj + cc[j]*hh[j+1];
359:   }

361:   /*
362:     compute the new plane rotation, and apply it to:
363:      1) the right-hand-side of the Hessenberg system (RS)
364:         note: it affects RS(it) and RS(it+1)
365:      2) the new column of the Hessenberg matrix
366:         note: it affects HH(it,it) which is currently pointed to
367:         by hh and HH(it+1, it) (*(hh+1))
368:     thus obtaining the updated value of the residual...
369:   */

371:   /* compute new plane rotation */

373:   if (!*hapend) {
374:     PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1])));
375:     if (delta == 0.0) {
376:       ksp->reason = KSP_DIVERGED_NULL;
377:       return(0);
378:     }

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

383:     hh[it]   = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1];
384:     rs[it+1] = -ss[it]*rs[it];
385:     rs[it]   = PetscConj(cc[it])*rs[it];
386:     *res     = PetscAbsScalar(rs[it+1]);
387:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
388:             another rotation matrix (so RH doesn't change).  The new residual is
389:             always the new sine term times the residual from last time (RS(it)),
390:             but now the new sine rotation would be zero...so the residual should
391:             be zero...so we will multiply "zero" by the last residual.  This might
392:             not be exactly what we want to do here -could just return "zero". */

394:     *res = 0.0;
395:   }
396:   return(0);
397: }

399: /*
400:    KSPBuildSolution_PGMRES

402:      Input Parameter:
403: .     ksp - the Krylov space object
404: .     ptr-

406:    Output Parameter:
407: .     result - the solution

409:    Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
410:    calls directly.

412: */
415: PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)
416: {
417:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)ksp->data;

421:   if (!ptr) {
422:     if (!pgmres->sol_temp) {
423:       VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);
424:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)pgmres->sol_temp);
425:     }
426:     ptr = pgmres->sol_temp;
427:   }
428:   if (!pgmres->nrs) {
429:     /* allocate the work area */
430:     PetscMalloc1(pgmres->max_k,&pgmres->nrs);
431:     PetscLogObjectMemory((PetscObject)ksp,pgmres->max_k*sizeof(PetscScalar));
432:   }

434:   KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);
435:   if (result) *result = ptr;
436:   return(0);
437: }

441: PetscErrorCode KSPSetFromOptions_PGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
442: {

446:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
447:   PetscOptionsHead(PetscOptionsObject,"KSP pipelined GMRES Options");
448:   PetscOptionsTail();
449:   return(0);
450: }

454: PetscErrorCode KSPReset_PGMRES(KSP ksp)
455: {

459:   KSPReset_GMRES(ksp);
460:   return(0);
461: }

463: /*MC
464:      KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.

466:    Options Database Keys:
467: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
468: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
469: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
470:                              vectors are allocated as needed)
471: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
472: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
473: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
474:                                    stability of the classical Gram-Schmidt  orthogonalization.
475: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

477:    Level: beginner

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

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

486:    Developer Notes: This object is subclassed off of KSPGMRES

488: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES, KSPPIPECG, KSPPIPECR,
489:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
490:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
491:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov()
492: M*/

496: PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp)
497: {
498:   KSP_PGMRES     *pgmres;

502:   PetscNewLog(ksp,&pgmres);

504:   ksp->data                              = (void*)pgmres;
505:   ksp->ops->buildsolution                = KSPBuildSolution_PGMRES;
506:   ksp->ops->setup                        = KSPSetUp_PGMRES;
507:   ksp->ops->solve                        = KSPSolve_PGMRES;
508:   ksp->ops->reset                        = KSPReset_PGMRES;
509:   ksp->ops->destroy                      = KSPDestroy_PGMRES;
510:   ksp->ops->view                         = KSPView_GMRES;
511:   ksp->ops->setfromoptions               = KSPSetFromOptions_PGMRES;
512:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
513:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

515:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
516:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);

518:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
519:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
520:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
521:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
522:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
523:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
524:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

526:   pgmres->nextra_vecs    = 1;
527:   pgmres->haptol         = 1.0e-30;
528:   pgmres->q_preallocate  = 0;
529:   pgmres->delta_allocate = PGMRES_DELTA_DIRECTIONS;
530:   pgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
531:   pgmres->nrs            = 0;
532:   pgmres->sol_temp       = 0;
533:   pgmres->max_k          = PGMRES_DEFAULT_MAXK;
534:   pgmres->Rsvd           = 0;
535:   pgmres->orthogwork     = 0;
536:   pgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
537:   return(0);
538: }