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