xref: /petsc/src/snes/utils/convest.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 #include <petscconvest.h>            /*I "petscconvest.h" I*/
2 #include <petscdmplex.h>
3 #include <petscds.h>
4 
5 #include <petsc/private/petscconvestimpl.h>
6 
7 static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
8 {
9   PetscInt c;
10   for (c = 0; c < Nc; ++c) u[c] = 0.0;
11   return 0;
12 }
13 
14 /*@
15   PetscConvEstDestroy - Destroys a PetscConvEst object
16 
17   Collective on PetscConvEst
18 
19   Input Parameter:
20 . ce - The PetscConvEst object
21 
22   Level: beginner
23 
24 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
25 @*/
26 PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
27 {
28   PetscErrorCode ierr;
29 
30   PetscFunctionBegin;
31   if (!*ce) PetscFunctionReturn(0);
32   PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1);
33   if (--((PetscObject)(*ce))->refct > 0) {
34     *ce = NULL;
35     PetscFunctionReturn(0);
36   }
37   ierr = PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);CHKERRQ(ierr);
38   ierr = PetscFree((*ce)->errors);CHKERRQ(ierr);
39   ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr);
40   PetscFunctionReturn(0);
41 }
42 
43 /*@
44   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
45 
46   Collective on PetscConvEst
47 
48   Input Parameters:
49 . ce - The PetscConvEst object
50 
51   Level: beginner
52 
53 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
54 @*/
55 PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
56 {
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr);
61   ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr);
62   ierr = PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL);CHKERRQ(ierr);
63   ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr);
64   ierr = PetscOptionsEnd();
65   PetscFunctionReturn(0);
66 }
67 
68 /*@
69   PetscConvEstView - Views a PetscConvEst object
70 
71   Collective on PetscConvEst
72 
73   Input Parameters:
74 + ce     - The PetscConvEst object
75 - viewer - The PetscViewer object
76 
77   Level: beginner
78 
79 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
80 @*/
81 PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
82 {
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr);
87   ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr);
88   PetscFunctionReturn(0);
89 }
90 
91 /*@
92   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
93 
94   Not collective
95 
96   Input Parameter:
97 . ce     - The PetscConvEst object
98 
99   Output Parameter:
100 . solver - The solver
101 
102   Level: intermediate
103 
104 .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
105 @*/
106 PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
107 {
108   PetscFunctionBegin;
109   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
110   PetscValidPointer(solver, 2);
111   *solver = ce->solver;
112   PetscFunctionReturn(0);
113 }
114 
115 /*@
116   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
117 
118   Not collective
119 
120   Input Parameters:
121 + ce     - The PetscConvEst object
122 - solver - The solver
123 
124   Level: intermediate
125 
126   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
127 
128 .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate()
129 @*/
130 PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
131 {
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
136   PetscValidHeader(solver, 2);
137   ce->solver = solver;
138   ierr = (*ce->ops->setsolver)(ce, solver);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 /*@
143   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
144 
145   Collective on PetscConvEst
146 
147   Input Parameters:
148 . ce - The PetscConvEst object
149 
150   Level: beginner
151 
152 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
153 @*/
154 PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
155 {
156   PetscDS        ds;
157   PetscInt       Nf, f;
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   ierr = DMGetDS(ce->idm, &ds);CHKERRQ(ierr);
162   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
163   ce->Nf = PetscMax(Nf, 1);
164   ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
165   ierr = PetscMalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr);
166   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
167   for (f = 0; f < Nf; ++f) {
168     ierr = PetscDSGetExactSolution(ds, f, &ce->exactSol[f], &ce->ctxs[f]);CHKERRQ(ierr);
169     if (!ce->exactSol[f]) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %D", f);
170   }
171   PetscFunctionReturn(0);
172 }
173 
174 PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
175 {
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
180   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
181   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
182   ierr = (*ce->ops->initguess)(ce, r, dm, u);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
192   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
193   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
194   PetscValidRealPointer(errors, 5);
195   ierr = (*ce->ops->computeerror)(ce, r, dm, u, errors);CHKERRQ(ierr);
196   PetscFunctionReturn(0);
197 }
198 
199 static PetscErrorCode PetscConvEstMonitor_Private(PetscConvEst ce, PetscInt r)
200 {
201   MPI_Comm       comm;
202   PetscInt       f;
203   PetscErrorCode ierr;
204 
205   PetscFunctionBegin;
206   if (ce->monitor) {
207     PetscReal *errors = &ce->errors[r*ce->Nf];
208 
209     ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
210     ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
211     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
212     for (f = 0; f < ce->Nf; ++f) {
213       if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
214       if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
215       else                     {ierr = PetscPrintf(comm, "%g", (double) errors[f]);CHKERRQ(ierr);}
216     }
217     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
218     ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
224 {
225   PetscClassId   id;
226   PetscErrorCode ierr;
227 
228   PetscFunctionBegin;
229   ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr);
230   if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
231   ierr = SNESGetDM((SNES) ce->solver, &ce->idm);CHKERRQ(ierr);
232   PetscFunctionReturn(0);
233 }
234 
235 static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
236 {
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
245 {
246   PetscErrorCode ierr;
247 
248   PetscFunctionBegin;
249   ierr = DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
254 {
255   SNES           snes = (SNES) ce->solver;
256   DM            *dm;
257   PetscObject    disc;
258   PetscReal     *x, *y, slope, intercept;
259   PetscInt      *dof, Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
260   void          *ctx;
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
265   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
266   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
267   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
268   ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr);
269   ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr);
270   /* Loop over meshes */
271   dm[0] = ce->idm;
272   for (r = 0; r <= Nr; ++r) {
273     Vec           u;
274     PetscLogStage stage;
275     char          stageName[PETSC_MAX_PATH_LEN];
276     const char   *dmname, *uname;
277 
278     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
279     ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);
280     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
281     if (r > 0) {
282       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
283       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
284       ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr);
285       ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr);
286       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
287       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
288       for (f = 0; f <= ce->Nf; ++f) {
289         PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
290 
291         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
292         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
293       }
294     }
295     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
296     /* Create solution */
297     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
298     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
299     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
300     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
301     /* Setup solver */
302     ierr = SNESReset(snes);CHKERRQ(ierr);
303     ierr = SNESSetDM(snes, dm[r]);CHKERRQ(ierr);
304     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
305     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
306     /* Create initial guess */
307     ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr);
308     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
309     ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
310     ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
311     ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
312     for (f = 0; f < ce->Nf; ++f) {
313       PetscSection s, fs;
314       PetscInt     lsize;
315 
316       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
317       ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr);
318       ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
319       ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr);
320       ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));CHKERRQ(ierr);
321       ierr = PetscLogEventSetDof(ce->event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr);
322       ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr);
323     }
324     /* Monitor */
325     ierr = PetscConvEstMonitor_Private(ce, r);CHKERRQ(ierr);
326     if (!r) {
327       /* PCReset() does not wipe out the level structure */
328       KSP ksp;
329       PC  pc;
330 
331       ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
332       ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
333       ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr);
334     }
335     /* Cleanup */
336     ierr = VecDestroy(&u);CHKERRQ(ierr);
337     ierr = PetscLogStagePop();CHKERRQ(ierr);
338   }
339   for (r = 1; r <= Nr; ++r) {
340     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
341   }
342   /* Fit convergence rate */
343   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
344   for (f = 0; f < ce->Nf; ++f) {
345     for (r = 0; r <= Nr; ++r) {
346       x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
347       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
348     }
349     ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
350     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
351     alpha[f] = -slope * dim;
352   }
353   ierr = PetscFree2(x, y);CHKERRQ(ierr);
354   ierr = PetscFree2(dm, dof);CHKERRQ(ierr);
355   /* Restore solver */
356   ierr = SNESReset(snes);CHKERRQ(ierr);
357   {
358     /* PCReset() does not wipe out the level structure */
359     KSP ksp;
360     PC  pc;
361 
362     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
363     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
364     ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr);
365     ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */
366   }
367   ierr = SNESSetDM(snes, ce->idm);CHKERRQ(ierr);
368   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
369   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373 /*@
374   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
375 
376   Not collective
377 
378   Input Parameter:
379 . ce   - The PetscConvEst object
380 
381   Output Parameter:
382 . alpha - The convergence rate for each field
383 
384   Note: The convergence rate alpha is defined by
385 $ || u_\Delta - u_exact || < C \Delta^alpha
386 where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
387 spatial resolution and \Delta t for the temporal resolution.
388 
389 We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
390 based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
391 
392   Options database keys:
393 + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
394 - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
395 
396   Level: intermediate
397 
398 .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
399 @*/
400 PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
401 {
402   PetscInt       f;
403   PetscErrorCode ierr;
404 
405   PetscFunctionBegin;
406   if (ce->event < 0) {ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);CHKERRQ(ierr);}
407   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
408   ierr = (*ce->ops->getconvrate)(ce, alpha);CHKERRQ(ierr);
409   PetscFunctionReturn(0);
410 }
411 
412 /*@
413   PetscConvEstRateView - Displays the convergence rate to a viewer
414 
415    Collective on SNES
416 
417    Parameter:
418 +  snes - iterative context obtained from SNESCreate()
419 .  alpha - the convergence rate for each field
420 -  viewer - the viewer to display the reason
421 
422    Options Database Keys:
423 .  -snes_convergence_estimate - print the convergence rate
424 
425    Level: developer
426 
427 .seealso: PetscConvEstGetRate()
428 @*/
429 PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
430 {
431   PetscBool      isAscii;
432   PetscErrorCode ierr;
433 
434   PetscFunctionBegin;
435   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
436   if (isAscii) {
437     PetscInt Nf = ce->Nf, f;
438 
439     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
440     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
441     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
442     for (f = 0; f < Nf; ++f) {
443       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
444       ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr);
445     }
446     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
447     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
448     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
449   }
450   PetscFunctionReturn(0);
451 }
452 
453 /*@
454   PetscConvEstCreate - Create a PetscConvEst object
455 
456   Collective
457 
458   Input Parameter:
459 . comm - The communicator for the PetscConvEst object
460 
461   Output Parameter:
462 . ce   - The PetscConvEst object
463 
464   Level: beginner
465 
466 .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
467 @*/
468 PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
469 {
470   PetscErrorCode ierr;
471 
472   PetscFunctionBegin;
473   PetscValidPointer(ce, 2);
474   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
475   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
476   (*ce)->monitor = PETSC_FALSE;
477   (*ce)->r       = 2.0;
478   (*ce)->Nr      = 4;
479   (*ce)->event   = -1;
480   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
481   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
482   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
483   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
484   PetscFunctionReturn(0);
485 }
486