xref: /petsc/src/ksp/ksp/impls/gmres/pgmres/pgmres.c (revision f2cacb80554376ba60c1c183c00cd40ea8a0a80d)
1 
2 /*
3     This file implements PGMRES (a Pipelined Generalized Minimal Residual method)
4 */
5 
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
9 
10 static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP,PetscInt,PetscBool*,PetscReal*);
11 static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
12 
13 /*
14 
15     KSPSetUp_PGMRES - Sets up the workspace needed by pgmres.
16 
17     This is called once, usually automatically by KSPSolve() or KSPSetUp(),
18     but can be called directly by KSPSetUp().
19 
20 */
21 static PetscErrorCode KSPSetUp_PGMRES(KSP ksp)
22 {
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = KSPSetUp_GMRES(ksp);CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 /*
31 
32     KSPPGMRESCycle - Run pgmres, possibly with restart.  Return residual
33                   history if requested.
34 
35     input parameters:
36 .        pgmres  - structure containing parameters and work areas
37 
38     output parameters:
39 .        itcount - number of iterations used.  If null, ignored.
40 .        converged - 0 if not converged
41 
42     Notes:
43     On entry, the value in vector VEC_VV(0) should be
44     the initial residual.
45 
46 
47  */
48 static PetscErrorCode KSPPGMRESCycle(PetscInt *itcount,KSP ksp)
49 {
50   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);
51   PetscReal      res_norm,res,newnorm;
52   PetscErrorCode ierr;
53   PetscInt       it     = 0,j,k;
54   PetscBool      hapend = PETSC_FALSE;
55 
56   PetscFunctionBegin;
57   if (itcount) *itcount = 0;
58   ierr   = VecNormalize(VEC_VV(0),&res_norm);CHKERRQ(ierr);
59   KSPCheckNorm(ksp,res_norm);
60   res    = res_norm;
61   *RS(0) = res_norm;
62 
63   /* check for the convergence */
64   ierr       = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
65   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
66   else ksp->rnorm = 0;
67   ierr       = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
68   pgmres->it = it-2;
69   ierr = KSPLogResidualHistory(ksp,ksp->rnorm);CHKERRQ(ierr);
70   ierr = KSPMonitor(ksp,ksp->its,ksp->rnorm);CHKERRQ(ierr);
71   if (!res) {
72     ksp->reason = KSP_CONVERGED_ATOL;
73     ierr        = PetscInfo(ksp,"Converged due to zero residual norm on entry\n");CHKERRQ(ierr);
74     PetscFunctionReturn(0);
75   }
76 
77   ierr = (*ksp->converged)(ksp,ksp->its,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
78   for (; !ksp->reason; it++) {
79     Vec Zcur,Znext;
80     if (pgmres->vv_allocated <= it + VEC_OFFSET + 1) {
81       ierr = KSPGMRESGetNewVectors(ksp,it+1);CHKERRQ(ierr);
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. */
87 
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       ierr = KSP_PCApplyBAorAB(ksp,Zcur,Znext,VEC_TEMP_MATOP);CHKERRQ(ierr);
90     }
91 
92     if (it > 1) {               /* Complete the pending reduction */
93       ierr           = VecNormEnd(VEC_VV(it-1),NORM_2,&newnorm);CHKERRQ(ierr);
94       *HH(it-1,it-2) = newnorm;
95     }
96     if (it > 0) {               /* Finish the reduction computing the latest column of H */
97       ierr = VecMDotEnd(Zcur,it,&(VEC_VV(0)),HH(0,it-1));CHKERRQ(ierr);
98     }
99 
100     if (it > 1) {
101       /* normalize the base vector from two iterations ago, basis is complete up to here */
102       ierr = VecScale(VEC_VV(it-1),1./ *HH(it-1,it-2));CHKERRQ(ierr);
103 
104       ierr       = KSPPGMRESUpdateHessenberg(ksp,it-2,&hapend,&res);CHKERRQ(ierr);
105       pgmres->it = it-2;
106       ksp->its++;
107       if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
108       else ksp->rnorm = 0;
109 
110       ierr = (*ksp->converged)(ksp,ksp->its,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
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         ierr = KSPLogResidualHistory(ksp,ksp->rnorm);CHKERRQ(ierr);
113         ierr = KSPMonitor(ksp,ksp->its,ksp->rnorm);CHKERRQ(ierr);
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       }
124 
125       if (!(it < pgmres->max_k+1 && ksp->its < ksp->max_it)) break;
126 
127       /* The it-2 column of H was not scaled when we computed Zcur, apply correction */
128       ierr = VecScale(Zcur,1./ *HH(it-1,it-2));CHKERRQ(ierr);
129       /* And Znext computed in this iteration was computed using the under-scaled Zcur */
130       ierr = VecScale(Znext,1./ *HH(it-1,it-2));CHKERRQ(ierr);
131 
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     }
138 
139     if (it > 0) {
140       PetscScalar *work;
141       if (!pgmres->orthogwork) {ierr = PetscMalloc1(pgmres->max_k + 2,&pgmres->orthogwork);CHKERRQ(ierr);}
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       ierr = VecMAXPY(Znext,it+1,work,&VEC_VV(0));CHKERRQ(ierr);
167       ierr = VecAXPY(Znext,-*HH(it-1,it-1),Zcur);CHKERRQ(ierr);
168 
169       /* Orthogonalize Zcur against existing basis vectors. */
170       for (k=0; k<it; k++) work[k] = -*HH(k,it-1);
171       ierr = VecMAXPY(Zcur,it,work,&VEC_VV(0));CHKERRQ(ierr);
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       ierr = VecNormBegin(VEC_VV(it),NORM_2,&newnorm);CHKERRQ(ierr);
175     }
176 
177     /* Compute column of H (to the diagonal, but not the subdiagonal) to be able to orthogonalize the newest vector. */
178     ierr = VecMDotBegin(Znext,it+1,&VEC_VV(0),HH(0,it));CHKERRQ(ierr);
179 
180     /* Start an asynchronous split-mode reduction, the result of the MDot and Norm will be collected on the next iteration. */
181     ierr = PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Znext));CHKERRQ(ierr);
182   }
183 
184   if (itcount) *itcount = it-1; /* Number of iterations actually completed. */
185 
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   ierr = KSPPGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-2);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 /*
197     KSPSolve_PGMRES - This routine applies the PGMRES method.
198 
199 
200    Input Parameter:
201 .     ksp - the Krylov space object that was set to use pgmres
202 
203    Output Parameter:
204 .     outits - number of iterations used
205 
206 */
207 static PetscErrorCode KSPSolve_PGMRES(KSP ksp)
208 {
209   PetscErrorCode ierr;
210   PetscInt       its,itcount;
211   KSP_PGMRES     *pgmres    = (KSP_PGMRES*)ksp->data;
212   PetscBool      guess_zero = ksp->guess_zero;
213 
214   PetscFunctionBegin;
215   if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
216   ierr     = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
217   ksp->its = 0;
218   ierr     = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
219 
220   itcount     = 0;
221   ksp->reason = KSP_CONVERGED_ITERATING;
222   while (!ksp->reason) {
223     ierr     = KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);CHKERRQ(ierr);
224     ierr     = KSPPGMRESCycle(&its,ksp);CHKERRQ(ierr);
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   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode KSPDestroy_PGMRES(KSP ksp)
237 {
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   ierr = KSPDestroy_GMRES(ksp);CHKERRQ(ierr);
242   PetscFunctionReturn(0);
243 }
244 
245 /*
246     KSPPGMRESBuildSoln - create the solution from the starting vector and the
247                       current iterates.
248 
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)
255 
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;
261   PetscErrorCode ierr;
262   PetscInt       k,j;
263   KSP_PGMRES     *pgmres = (KSP_PGMRES*)(ksp->data);
264 
265   PetscFunctionBegin;
266   /* Solve for solution vector that minimizes the residual */
267 
268   if (it < 0) {                                 /* no pgmres steps have been performed */
269     ierr = VecCopy(vguess,vdest);CHKERRQ(ierr); /* VecCopy() is smart, exits immediately if vguess == vdest */
270     PetscFunctionReturn(0);
271   }
272 
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;
277 
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   }
283 
284   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
285   ierr = VecZeroEntries(VEC_TEMP);CHKERRQ(ierr);
286   ierr = VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));CHKERRQ(ierr);
287   ierr = KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);CHKERRQ(ierr);
288   /* add solution to previous solution */
289   if (vdest == vguess) {
290     ierr = VecAXPY(vdest,1.0,VEC_TEMP);CHKERRQ(ierr);
291   } else {
292     ierr = VecWAXPY(vdest,1.0,VEC_TEMP,vguess);CHKERRQ(ierr);
293   }
294   PetscFunctionReturn(0);
295 }
296 
297 /*
298 
299     KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
300                             Return new residual.
301 
302     input parameters:
303 
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.
308 
309     output parameters:
310 .        res - the new residual
311 
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);
322   PetscErrorCode ierr;
323 
324   PetscFunctionBegin;
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 */
329 
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];
332 
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     ierr    = 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)));CHKERRQ(ierr);
337     *hapend = PETSC_TRUE;
338   }
339 
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!) */
344 
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   }
350 
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   */
360 
361   /* compute new plane rotation */
362 
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       PetscFunctionReturn(0);
368     }
369 
370     cc[it] = hh[it] / delta;    /* new cosine value */
371     ss[it] = hh[it+1] / delta;  /* new sine value */
372 
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". */
383 
384     *res = 0.0;
385   }
386   PetscFunctionReturn(0);
387 }
388 
389 /*
390    KSPBuildSolution_PGMRES
391 
392      Input Parameter:
393 .     ksp - the Krylov space object
394 .     ptr-
395 
396    Output Parameter:
397 .     result - the solution
398 
399    Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle
400    calls directly.
401 
402 */
403 PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result)
404 {
405   KSP_PGMRES     *pgmres = (KSP_PGMRES*)ksp->data;
406   PetscErrorCode ierr;
407 
408   PetscFunctionBegin;
409   if (!ptr) {
410     if (!pgmres->sol_temp) {
411       ierr = VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);CHKERRQ(ierr);
412       ierr = PetscLogObjectParent((PetscObject)ksp,(PetscObject)pgmres->sol_temp);CHKERRQ(ierr);
413     }
414     ptr = pgmres->sol_temp;
415   }
416   if (!pgmres->nrs) {
417     /* allocate the work area */
418     ierr = PetscMalloc1(pgmres->max_k,&pgmres->nrs);CHKERRQ(ierr);
419     ierr = PetscLogObjectMemory((PetscObject)ksp,pgmres->max_k*sizeof(PetscScalar));CHKERRQ(ierr);
420   }
421 
422   ierr = KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);CHKERRQ(ierr);
423   if (result) *result = ptr;
424   PetscFunctionReturn(0);
425 }
426 
427 PetscErrorCode KSPSetFromOptions_PGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
428 {
429   PetscErrorCode ierr;
430 
431   PetscFunctionBegin;
432   ierr = KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);CHKERRQ(ierr);
433   ierr = PetscOptionsHead(PetscOptionsObject,"KSP pipelined GMRES Options");CHKERRQ(ierr);
434   ierr = PetscOptionsTail();CHKERRQ(ierr);
435   PetscFunctionReturn(0);
436 }
437 
438 PetscErrorCode KSPReset_PGMRES(KSP ksp)
439 {
440   PetscErrorCode ierr;
441 
442   PetscFunctionBegin;
443   ierr = KSPReset_GMRES(ksp);CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 /*MC
448      KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method.
449 
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
460 
461    Level: beginner
462 
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.
466 
467    Reference:
468    Ghysels, Ashby, Meerbergen, Vanroose, Hiding global communication latencies in the GMRES algorithm on massively parallel machines, 2012.
469 
470    Developer Notes:
471     This object is subclassed off of KSPGMRES
472 
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*/
478 
479 PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp)
480 {
481   KSP_PGMRES     *pgmres;
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   ierr = PetscNewLog(ksp,&pgmres);CHKERRQ(ierr);
486 
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;
497 
498   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
499   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);CHKERRQ(ierr);
500   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);CHKERRQ(ierr);
501 
502   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);CHKERRQ(ierr);
503   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);CHKERRQ(ierr);
504   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);CHKERRQ(ierr);
505   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);CHKERRQ(ierr);
506   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);CHKERRQ(ierr);
507   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);CHKERRQ(ierr);
508   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);CHKERRQ(ierr);
509 
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   PetscFunctionReturn(0);
522 }
523