Actual source code: pgmres.c

petsc-3.3-p5 2012-12-01
  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:   res    = res_norm;
 64:   *RS(0) = res_norm;

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

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

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

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

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

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

111:       (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
112:       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. */
113:         KSPLogResidualHistory(ksp,res);
114:         KSPMonitor(ksp,ksp->its,res);
115:       }
116:       if (ksp->reason) break;
117:       /* Catch error in happy breakdown and signal convergence and break from loop */
118:       if (hapend) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down, but convergence was not indicated. Residual norm = %G",res);
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) {PetscMalloc((pgmres->max_k + 2)*sizeof(PetscScalar),&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(((PetscObject)Znext)->comm);
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.


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

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

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

211:   if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
212:   PetscObjectTakeAccess(ksp);
213:   ksp->its = 0;
214:   PetscObjectGrantAccess(ksp);

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

234: static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)
235: {

239:   KSPDestroy_GMRES(ksp);
240:   return(0);
241: }

243: /*
244:     KSPPGMRESBuildSoln - create the solution from the starting vector and the
245:                       current iterates.

247:     Input parameters:
248:         nrs - work area of size it + 1.
249:         vguess  - index of initial guess
250:         vdest - index of result.  Note that vguess may == vdest (replace
251:                 guess with the solution).
252:         it - HH upper triangular part is a block of size (it+1) x (it+1)

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

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

268:   if (it < 0) {                                 /* no pgmres steps have been performed */
269:     VecCopy(vguess,vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
270:     return(0);
271:   }

273:   /* solve the upper triangular system - RS is the right side and HH is
274:      the upper triangular matrix  - put soln in nrs */
275:   if (*HH(it,it) != 0.0) {
276:     nrs[it] = *RS(it) / *HH(it,it);
277:   } else {
278:     nrs[it] = 0.0;
279:   }
280:   for (k=it-1; k>=0; k--) {
281:     tt  = *RS(k);
282:     for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
283:     nrs[k]   = tt / *HH(k,k);
284:   }

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

299: /*

301:     KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
302:                             Return new residual.

304:     input parameters:

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

311:     output parameters:
312: .        res - the new residual

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

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

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

337:   /* check for the happy breakdown */
338:   hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pgmres->haptol);
339:   if (PetscAbsScalar(hh[it+1]) < hapbnd) {
340:     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)));
341:     *hapend = PETSC_TRUE;
342:   }

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

349:   for (j=0; j<it; j++) {
350:     PetscScalar hhj = hh[j];
351:     hh[j]   = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
352:     hh[j+1] =          -ss[j] *hhj + cc[j]*hh[j+1];
353:   }

355:   /*
356:     compute the new plane rotation, and apply it to:
357:      1) the right-hand-side of the Hessenberg system (RS)
358:         note: it affects RS(it) and RS(it+1)
359:      2) the new column of the Hessenberg matrix
360:         note: it affects HH(it,it) which is currently pointed to
361:         by hh and HH(it+1, it) (*(hh+1))
362:     thus obtaining the updated value of the residual...
363:   */

365:   /* compute new plane rotation */

367:   if (!*hapend) {
368:     PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1])));
369:     if (delta == 0.0) {
370:       ksp->reason = KSP_DIVERGED_NULL;
371:       return(0);
372:     }

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

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

388:     *res = 0.0;
389:   }
390:   return(0);
391: }

393: /*
394:    KSPBuildSolution_PGMRES

396:      Input Parameter:
397: .     ksp - the Krylov space object
398: .     ptr-

400:    Output Parameter:
401: .     result - the solution

403:    Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
404:    calls directly.

406: */
409: PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)
410: {
411:   KSP_PGMRES     *pgmres = (KSP_PGMRES *)ksp->data;

415:   if (!ptr) {
416:     if (!pgmres->sol_temp) {
417:       VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);
418:       PetscLogObjectParent(ksp,pgmres->sol_temp);
419:     }
420:     ptr = pgmres->sol_temp;
421:   }
422:   if (!pgmres->nrs) {
423:     /* allocate the work area */
424:     PetscMalloc(pgmres->max_k*sizeof(PetscScalar),&pgmres->nrs);
425:     PetscLogObjectMemory(ksp,pgmres->max_k*sizeof(PetscScalar));
426:   }

428:   KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);
429:   if (result) *result = ptr;

431:   return(0);
432: }

436: PetscErrorCode KSPSetFromOptions_PGMRES(KSP ksp)
437: {

441:   KSPSetFromOptions_GMRES(ksp);
442:   PetscOptionsHead("KSP pipelined GMRES Options");
443:   PetscOptionsTail();
444:   return(0);
445: }

449: PetscErrorCode KSPReset_PGMRES(KSP ksp)
450: {

454:   KSPReset_GMRES(ksp);
455:   return(0);
456: }

458: /*MC
459:      KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.

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

472:    Level: beginner

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

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

479: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
480:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
481:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
482:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov()
483: M*/

488: {
489:   KSP_PGMRES     *pgmres;

493:   PetscNewLog(ksp,KSP_PGMRES,&pgmres);
494:   ksp->data                              = (void*)pgmres;
495:   ksp->ops->buildsolution                = KSPBuildSolution_PGMRES;
496:   ksp->ops->setup                        = KSPSetUp_PGMRES;
497:   ksp->ops->solve                        = KSPSolve_PGMRES;
498:   ksp->ops->reset                        = KSPReset_PGMRES;
499:   ksp->ops->destroy                      = KSPDestroy_PGMRES;
500:   ksp->ops->view                         = KSPView_GMRES;
501:   ksp->ops->setfromoptions               = KSPSetFromOptions_PGMRES;
502:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
503:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

505:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
506:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);

508:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
509:                                     "KSPGMRESSetPreAllocateVectors_GMRES",
510:                                      KSPGMRESSetPreAllocateVectors_GMRES);
511:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
512:                                     "KSPGMRESSetOrthogonalization_GMRES",
513:                                      KSPGMRESSetOrthogonalization_GMRES);
514:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
515:                                     "KSPGMRESGetOrthogonalization_GMRES",
516:                                      KSPGMRESGetOrthogonalization_GMRES);
517:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
518:                                     "KSPGMRESSetRestart_GMRES",
519:                                      KSPGMRESSetRestart_GMRES);
520:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
521:                                     "KSPGMRESGetRestart_GMRES",
522:                                      KSPGMRESGetRestart_GMRES);
523:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
524:                                     "KSPGMRESSetCGSRefinementType_GMRES",
525:                                      KSPGMRESSetCGSRefinementType_GMRES);
526:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
527:                                     "KSPGMRESGetCGSRefinementType_GMRES",
528:                                      KSPGMRESGetCGSRefinementType_GMRES);

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