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: PetscObjectAMSTakeAccess((PetscObject)ksp);
68: ksp->rnorm = res;
69: PetscObjectAMSGrantAccess((PetscObject)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) {
119: 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",res);
120: else {
121: ksp->reason = KSP_DIVERGED_BREAKDOWN;
122: break;
123: }
124: }
126: if (!(it < pgmres->max_k+1 && ksp->its < ksp->max_it)) break;
128: /* The it-2 column of H was not scaled when we computed Zcur, apply correction */
129: VecScale(Zcur,1./ *HH(it-1,it-2));
130: /* And Znext computed in this iteration was computed using the under-scaled Zcur */
131: VecScale(Znext,1./ *HH(it-1,it-2));
133: /* 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. */
134: for (k=0; k<it; k++) *HH(k,it-1) /= *HH(it-1,it-2);
135: /* When Zcur was projected against the Krylov basis, VV(it-1) was still not normalized, so fix that too. This
136: * column is complete except for HH(it,it-1) which we won't know until the next iteration. */
137: *HH(it-1,it-1) /= *HH(it-1,it-2);
138: }
140: if (it > 0) {
141: PetscScalar *work;
142: if (!pgmres->orthogwork) {PetscMalloc((pgmres->max_k + 2)*sizeof(PetscScalar),&pgmres->orthogwork);}
143: work = pgmres->orthogwork;
144: /* Apply correction computed by the VecMDot in the last iteration to Znext. The original form is
145: *
146: * Znext -= sum_{j=0}^{i-1} Z[j+1] * H[j,i-1]
147: *
148: * where
149: *
150: * Z[j] = sum_{k=0}^j V[k] * H[k,j-1]
151: *
152: * substituting
153: *
154: * Znext -= sum_{j=0}^{i-1} sum_{k=0}^{j+1} V[k] * H[k,j] * H[j,i-1]
155: *
156: * rearranging the iteration space from row-column to column-row
157: *
158: * Znext -= sum_{k=0}^i sum_{j=k-1}^{i-1} V[k] * H[k,j] * H[j,i-1]
159: *
160: * Note that column it-1 of HH is correct. For all previous columns, we must look at HES because HH has already
161: * been transformed to upper triangular form.
162: */
163: for (k=0; k<it+1; k++) {
164: work[k] = 0;
165: for (j=PetscMax(0,k-1); j<it-1; j++) work[k] -= *HES(k,j) * *HH(j,it-1);
166: }
167: VecMAXPY(Znext,it+1,work,&VEC_VV(0));
168: VecAXPY(Znext,-*HH(it-1,it-1),Zcur);
170: /* Orthogonalize Zcur against existing basis vectors. */
171: for (k=0; k<it; k++) work[k] = -*HH(k,it-1);
172: VecMAXPY(Zcur,it,work,&VEC_VV(0));
173: /* Zcur is now orthogonal, and will be referred to as VEC_VV(it) again, though it is still not normalized. */
174: /* Begin computing the norm of the new vector, will be normalized after the MatMult in the next iteration. */
175: VecNormBegin(VEC_VV(it),NORM_2,&newnorm);
176: }
178: /* Compute column of H (to the diagonal, but not the subdiagonal) to be able to orthogonalize the newest vector. */
179: VecMDotBegin(Znext,it+1,&VEC_VV(0),HH(0,it));
181: /* Start an asynchronous split-mode reduction, the result of the MDot and Norm will be collected on the next iteration. */
182: PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Znext));
183: }
185: if (itcount) *itcount = it-1; /* Number of iterations actually completed. */
187: /*
188: Down here we have to solve for the "best" coefficients of the Krylov
189: columns, add the solution values together, and possibly unwind the
190: preconditioning from the solution
191: */
192: /* Form the solution (or the solution so far) */
193: KSPPGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-2);
194: return(0);
195: }
197: /*
198: KSPSolve_PGMRES - This routine applies the PGMRES method.
201: Input Parameter:
202: . ksp - the Krylov space object that was set to use pgmres
204: Output Parameter:
205: . outits - number of iterations used
207: */
210: static PetscErrorCode KSPSolve_PGMRES(KSP ksp)211: {
213: PetscInt its,itcount;
214: KSP_PGMRES *pgmres = (KSP_PGMRES*)ksp->data;
215: PetscBool guess_zero = ksp->guess_zero;
218: if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
219: PetscObjectAMSTakeAccess((PetscObject)ksp);
220: ksp->its = 0;
221: PetscObjectAMSGrantAccess((PetscObject)ksp);
223: itcount = 0;
224: ksp->reason = KSP_CONVERGED_ITERATING;
225: while (!ksp->reason) {
226: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
227: KSPPGMRESCycle(&its,ksp);
228: itcount += its;
229: if (itcount >= ksp->max_it) {
230: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
231: break;
232: }
233: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
234: }
235: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
236: return(0);
237: }
241: static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)242: {
246: KSPDestroy_GMRES(ksp);
247: return(0);
248: }
250: /*
251: KSPPGMRESBuildSoln - create the solution from the starting vector and the
252: current iterates.
254: Input parameters:
255: nrs - work area of size it + 1.
256: vguess - index of initial guess
257: vdest - index of result. Note that vguess may == vdest (replace
258: guess with the solution).
259: it - HH upper triangular part is a block of size (it+1) x (it+1)
261: This is an internal routine that knows about the PGMRES internals.
262: */
265: static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)266: {
267: PetscScalar tt;
269: PetscInt k,j;
270: KSP_PGMRES *pgmres = (KSP_PGMRES*)(ksp->data);
273: /* Solve for solution vector that minimizes the residual */
275: if (it < 0) { /* no pgmres steps have been performed */
276: VecCopy(vguess,vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
277: return(0);
278: }
280: /* solve the upper triangular system - RS is the right side and HH is
281: the upper triangular matrix - put soln in nrs */
282: if (*HH(it,it) != 0.0) nrs[it] = *RS(it) / *HH(it,it);
283: else nrs[it] = 0.0;
285: for (k=it-1; k>=0; k--) {
286: tt = *RS(k);
287: for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
288: nrs[k] = tt / *HH(k,k);
289: }
291: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
292: VecZeroEntries(VEC_TEMP);
293: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
294: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
295: /* add solution to previous solution */
296: if (vdest == vguess) {
297: VecAXPY(vdest,1.0,VEC_TEMP);
298: } else {
299: VecWAXPY(vdest,1.0,VEC_TEMP,vguess);
300: }
301: return(0);
302: }
304: /*
306: KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
307: Return new residual.
309: input parameters:
311: . ksp - Krylov space object
312: . it - plane rotations are applied to the (it+1)th column of the
313: modified hessenberg (i.e. HH(:,it))
314: . hapend - PETSC_FALSE not happy breakdown ending.
316: output parameters:
317: . res - the new residual
319: */
322: /*
323: . it - column of the Hessenberg that is complete, PGMRES is actually computing two columns ahead of this
324: */
325: static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res)326: {
327: PetscScalar *hh,*cc,*ss,*rs;
328: PetscInt j;
329: PetscReal hapbnd;
330: KSP_PGMRES *pgmres = (KSP_PGMRES*)(ksp->data);
334: hh = HH(0,it); /* pointer to beginning of column to update */
335: cc = CC(0); /* beginning of cosine rotations */
336: ss = SS(0); /* beginning of sine rotations */
337: rs = RS(0); /* right hand side of least squares system */
339: /* The Hessenberg matrix is now correct through column it, save that form for possible spectral analysis */
340: for (j=0; j<=it+1; j++) *HES(j,it) = hh[j];
342: /* check for the happy breakdown */
343: hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pgmres->haptol);
344: if (PetscAbsScalar(hh[it+1]) < hapbnd) {
345: 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)));
346: *hapend = PETSC_TRUE;
347: }
349: /* Apply all the previously computed plane rotations to the new column
350: of the Hessenberg matrix */
351: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta),
352: and some refs have [c s ; -conj(s) c] (don't be confused!) */
354: for (j=0; j<it; j++) {
355: PetscScalar hhj = hh[j];
356: hh[j] = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
357: hh[j+1] = -ss[j] *hhj + cc[j]*hh[j+1];
358: }
360: /*
361: compute the new plane rotation, and apply it to:
362: 1) the right-hand-side of the Hessenberg system (RS)
363: note: it affects RS(it) and RS(it+1)
364: 2) the new column of the Hessenberg matrix
365: note: it affects HH(it,it) which is currently pointed to
366: by hh and HH(it+1, it) (*(hh+1))
367: thus obtaining the updated value of the residual...
368: */
370: /* compute new plane rotation */
372: if (!*hapend) {
373: PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1])));
374: if (delta == 0.0) {
375: ksp->reason = KSP_DIVERGED_NULL;
376: return(0);
377: }
379: cc[it] = hh[it] / delta; /* new cosine value */
380: ss[it] = hh[it+1] / delta; /* new sine value */
382: hh[it] = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1];
383: rs[it+1] = -ss[it]*rs[it];
384: rs[it] = PetscConj(cc[it])*rs[it];
385: *res = PetscAbsScalar(rs[it+1]);
386: } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
387: another rotation matrix (so RH doesn't change). The new residual is
388: always the new sine term times the residual from last time (RS(it)),
389: but now the new sine rotation would be zero...so the residual should
390: be zero...so we will multiply "zero" by the last residual. This might
391: not be exactly what we want to do here -could just return "zero". */
393: *res = 0.0;
394: }
395: return(0);
396: }
398: /*
399: KSPBuildSolution_PGMRES
401: Input Parameter:
402: . ksp - the Krylov space object
403: . ptr-
405: Output Parameter:
406: . result - the solution
408: Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
409: calls directly.
411: */
414: PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)415: {
416: KSP_PGMRES *pgmres = (KSP_PGMRES*)ksp->data;
420: if (!ptr) {
421: if (!pgmres->sol_temp) {
422: VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);
423: PetscLogObjectParent(ksp,pgmres->sol_temp);
424: }
425: ptr = pgmres->sol_temp;
426: }
427: if (!pgmres->nrs) {
428: /* allocate the work area */
429: PetscMalloc(pgmres->max_k*sizeof(PetscScalar),&pgmres->nrs);
430: PetscLogObjectMemory(ksp,pgmres->max_k*sizeof(PetscScalar));
431: }
433: KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);
434: if (result) *result = ptr;
435: return(0);
436: }
440: PetscErrorCode KSPSetFromOptions_PGMRES(KSP ksp)441: {
445: KSPSetFromOptions_GMRES(ksp);
446: PetscOptionsHead("KSP pipelined GMRES Options");
447: PetscOptionsTail();
448: return(0);
449: }
453: PetscErrorCode KSPReset_PGMRES(KSP ksp)454: {
458: KSPReset_GMRES(ksp);
459: return(0);
460: }
462: /*MC
463: KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.
465: Options Database Keys:
466: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
467: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
468: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
469: vectors are allocated as needed)
470: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
471: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
472: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
473: stability of the classical Gram-Schmidt orthogonalization.
474: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
476: Level: beginner
478: Notes:
479: MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
480: See the FAQ on the PETSc website for details.
482: Reference:
483: Ghysels, Ashby, Meerbergen, Vanroose, Hiding global communication latencies in the GMRES algorithm on massively parallel machines, 2012.
485: Developer Notes: This object is subclassed off of KSPGMRES487: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES, KSPPIPECG, KSPPIPECR,
488: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
489: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
490: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov()
491: M*/
495: PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp)496: {
497: KSP_PGMRES *pgmres;
501: PetscNewLog(ksp,KSP_PGMRES,&pgmres);
503: ksp->data = (void*)pgmres;
504: ksp->ops->buildsolution = KSPBuildSolution_PGMRES;
505: ksp->ops->setup = KSPSetUp_PGMRES;
506: ksp->ops->solve = KSPSolve_PGMRES;
507: ksp->ops->reset = KSPReset_PGMRES;
508: ksp->ops->destroy = KSPDestroy_PGMRES;
509: ksp->ops->view = KSPView_GMRES;
510: ksp->ops->setfromoptions = KSPSetFromOptions_PGMRES;
511: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
512: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
514: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
515: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
517: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
518: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
519: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
520: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
521: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
522: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
523: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
525: pgmres->nextra_vecs = 1;
526: pgmres->haptol = 1.0e-30;
527: pgmres->q_preallocate = 0;
528: pgmres->delta_allocate = PGMRES_DELTA_DIRECTIONS;
529: pgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
530: pgmres->nrs = 0;
531: pgmres->sol_temp = 0;
532: pgmres->max_k = PGMRES_DEFAULT_MAXK;
533: pgmres->Rsvd = 0;
534: pgmres->orthogwork = 0;
535: pgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
536: return(0);
537: }