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