Actual source code: pgmres.c

petsc-3.14.6 2021-03-30
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>
  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: {

 26:   KSPSetUp_GMRES(ksp);
 27:   return(0);
 28: }

 30: /*

 32:     KSPPGMRESCycle - Run pgmres, possibly with restart.  Return residual
 33:                   history if requested.

 35:     input parameters:
 36: .        pgmres  - structure containing parameters and work areas

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

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


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

 57:   if (itcount) *itcount = 0;
 58:   VecNormalize(VEC_VV(0),&res_norm);
 59:   KSPCheckNorm(ksp,res_norm);
 60:   res    = res_norm;
 61:   *RS(0) = res_norm;

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

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

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

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

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

104:       KSPPGMRESUpdateHessenberg(ksp,it-2,&hapend,&res);
105:       pgmres->it = it-2;
106:       ksp->its++;
107:       if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
108:       else ksp->rnorm = 0;

110:       (*ksp->converged)(ksp,ksp->its,ksp->rnorm,&ksp->reason,ksp->cnvP);
111:       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. */
112:         KSPLogResidualHistory(ksp,ksp->rnorm);
113:         KSPMonitor(ksp,ksp->its,ksp->rnorm);
114:       }
115:       if (ksp->reason) break;
116:       /* Catch error in happy breakdown and signal convergence and break from loop */
117:       if (hapend) {
118:         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);
119:         else {
120:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
121:           break;
122:         }
123:       }

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

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

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

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

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

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

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

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

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

196: /*
197:     KSPSolve_PGMRES - This routine applies the PGMRES method.


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

203:    Output Parameter:
204: .     outits - number of iterations used

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

215:   if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
216:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
217:   ksp->its = 0;
218:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

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

236: static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)
237: {

241:   KSPDestroy_GMRES(ksp);
242:   return(0);
243: }

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

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

256:      This is an internal routine that knows about the PGMRES internals.
257:  */
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) nrs[it] = *RS(it) / *HH(it,it);
276:   else nrs[it] = 0.0;

278:   for (k=it-1; k>=0; k--) {
279:     tt = *RS(k);
280:     for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
281:     nrs[k] = tt / *HH(k,k);
282:   }

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

297: /*

299:     KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
300:                             Return new residual.

302:     input parameters:

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

309:     output parameters:
310: .        res - the new residual

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

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

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

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

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

345:   for (j=0; j<it; j++) {
346:     PetscScalar hhj = hh[j];
347:     hh[j]   = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
348:     hh[j+1] =          -ss[j] *hhj + cc[j]*hh[j+1];
349:   }

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

361:   /* compute new plane rotation */

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

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

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

384:     *res = 0.0;
385:   }
386:   return(0);
387: }

389: /*
390:    KSPBuildSolution_PGMRES

392:      Input Parameter:
393: .     ksp - the Krylov space object
394: .     ptr-

396:    Output Parameter:
397: .     result - the solution

399:    Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
400:    calls directly.

402: */
403: PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)
404: {
405:   KSP_PGMRES     *pgmres = (KSP_PGMRES*)ksp->data;

409:   if (!ptr) {
410:     if (!pgmres->sol_temp) {
411:       VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);
412:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)pgmres->sol_temp);
413:     }
414:     ptr = pgmres->sol_temp;
415:   }
416:   if (!pgmres->nrs) {
417:     /* allocate the work area */
418:     PetscMalloc1(pgmres->max_k,&pgmres->nrs);
419:     PetscLogObjectMemory((PetscObject)ksp,pgmres->max_k*sizeof(PetscScalar));
420:   }

422:   KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);
423:   if (result) *result = ptr;
424:   return(0);
425: }

427: PetscErrorCode KSPSetFromOptions_PGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
428: {

432:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
433:   PetscOptionsHead(PetscOptionsObject,"KSP pipelined GMRES Options");
434:   PetscOptionsTail();
435:   return(0);
436: }

438: PetscErrorCode KSPReset_PGMRES(KSP ksp)
439: {

443:   KSPReset_GMRES(ksp);
444:   return(0);
445: }

447: /*MC
448:      KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.

450:    Options Database Keys:
451: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
452: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
453: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
454:                              vectors are allocated as needed)
455: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
456: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
457: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
458:                                    stability of the classical Gram-Schmidt  orthogonalization.
459: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated

461:    Level: beginner

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

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

470:    Developer Notes:
471:     This object is subclassed off of KSPGMRES

473: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES, KSPPIPECG, KSPPIPECR,
474:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
475:            KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
476:            KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(),  KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov()
477: M*/

479: PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp)
480: {
481:   KSP_PGMRES     *pgmres;

485:   PetscNewLog(ksp,&pgmres);

487:   ksp->data                              = (void*)pgmres;
488:   ksp->ops->buildsolution                = KSPBuildSolution_PGMRES;
489:   ksp->ops->setup                        = KSPSetUp_PGMRES;
490:   ksp->ops->solve                        = KSPSolve_PGMRES;
491:   ksp->ops->reset                        = KSPReset_PGMRES;
492:   ksp->ops->destroy                      = KSPDestroy_PGMRES;
493:   ksp->ops->view                         = KSPView_GMRES;
494:   ksp->ops->setfromoptions               = KSPSetFromOptions_PGMRES;
495:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
496:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

498:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
499:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
500:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);

502:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
503:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
504:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
505:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
506:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
507:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
508:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

510:   pgmres->nextra_vecs    = 1;
511:   pgmres->haptol         = 1.0e-30;
512:   pgmres->q_preallocate  = 0;
513:   pgmres->delta_allocate = PGMRES_DELTA_DIRECTIONS;
514:   pgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
515:   pgmres->nrs            = NULL;
516:   pgmres->sol_temp       = NULL;
517:   pgmres->max_k          = PGMRES_DEFAULT_MAXK;
518:   pgmres->Rsvd           = NULL;
519:   pgmres->orthogwork     = NULL;
520:   pgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
521:   return(0);
522: }