xref: /petsc/src/snes/utils/convest.c (revision d5c9c0c4eebc2f2a01a1bd0c86fca87e2acd2a03)
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   PetscInt       Nf, f, Nds, s;
157   PetscErrorCode ierr;
158 
159   PetscFunctionBegin;
160   ierr = DMGetNumFields(ce->idm, &Nf);CHKERRQ(ierr);
161   ce->Nf = PetscMax(Nf, 1);
162   ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
163   ierr = PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr);
164   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
165   ierr = DMGetNumDS(ce->idm, &Nds);CHKERRQ(ierr);
166   for (s = 0; s < Nds; ++s) {
167     PetscDS         ds;
168     DMLabel         label;
169     IS              fieldIS;
170     const PetscInt *fields;
171     PetscInt        dsNf;
172 
173     ierr = DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
174     ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
175     if (fieldIS) {ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);}
176     for (f = 0; f < dsNf; ++f) {
177       const PetscInt field = fields[f];
178       ierr = PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]);CHKERRQ(ierr);
179     }
180     if (fieldIS) {ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);}
181   }
182   for (f = 0; f < Nf; ++f) {
183     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);
184   }
185   PetscFunctionReturn(0);
186 }
187 
188 PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
189 {
190   PetscErrorCode ierr;
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   ierr = (*ce->ops->initguess)(ce, r, dm, u);CHKERRQ(ierr);
197   PetscFunctionReturn(0);
198 }
199 
200 PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
206   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
207   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
208   PetscValidRealPointer(errors, 5);
209   ierr = (*ce->ops->computeerror)(ce, r, dm, u, errors);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode PetscConvEstMonitor_Private(PetscConvEst ce, PetscInt r)
214 {
215   MPI_Comm       comm;
216   PetscInt       f;
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   if (ce->monitor) {
221     PetscReal *errors = &ce->errors[r*ce->Nf];
222 
223     ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
224     ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
225     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
226     for (f = 0; f < ce->Nf; ++f) {
227       if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
228       if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
229       else                     {ierr = PetscPrintf(comm, "%g", (double) errors[f]);CHKERRQ(ierr);}
230     }
231     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
232     ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
233   }
234   PetscFunctionReturn(0);
235 }
236 
237 static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
238 {
239   PetscClassId   id;
240   PetscErrorCode ierr;
241 
242   PetscFunctionBegin;
243   ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr);
244   if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
245   ierr = SNESGetDM((SNES) ce->solver, &ce->idm);CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
250 {
251   PetscErrorCode ierr;
252 
253   PetscFunctionBegin;
254   ierr = DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr);
255   PetscFunctionReturn(0);
256 }
257 
258 static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
259 {
260   PetscErrorCode ierr;
261 
262   PetscFunctionBegin;
263   ierr = DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr);
264   PetscFunctionReturn(0);
265 }
266 
267 static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
268 {
269   SNES           snes = (SNES) ce->solver;
270   DM            *dm;
271   PetscObject    disc;
272   PetscReal     *x, *y, slope, intercept;
273   PetscInt      *dof, Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
274   void          *ctx;
275   PetscErrorCode ierr;
276 
277   PetscFunctionBegin;
278   if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
279   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
280   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
281   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
282   ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr);
283   ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr);
284   /* Loop over meshes */
285   dm[0] = ce->idm;
286   for (r = 0; r <= Nr; ++r) {
287     Vec           u;
288 #if defined(PETSC_USE_LOG)
289     PetscLogStage stage;
290 #endif
291     char          stageName[PETSC_MAX_PATH_LEN];
292     const char   *dmname, *uname;
293 
294     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
295     ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);
296     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
297     if (r > 0) {
298       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
299       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
300       ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr);
301       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
302       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
303       for (f = 0; f <= ce->Nf; ++f) {
304         PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
305 
306         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
307         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
308       }
309     }
310     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
311     /* Create solution */
312     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
313     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
314     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
315     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
316     /* Setup solver */
317     ierr = SNESReset(snes);CHKERRQ(ierr);
318     ierr = SNESSetDM(snes, dm[r]);CHKERRQ(ierr);
319     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
320     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
321     /* Create initial guess */
322     ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr);
323     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
324     ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
325     ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
326     ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
327     for (f = 0; f < ce->Nf; ++f) {
328       PetscSection s, fs;
329       PetscInt     lsize;
330 
331       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
332       ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr);
333       ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
334       ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr);
335       ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));CHKERRQ(ierr);
336       ierr = PetscLogEventSetDof(ce->event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr);
337       ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr);
338     }
339     /* Monitor */
340     ierr = PetscConvEstMonitor_Private(ce, r);CHKERRQ(ierr);
341     if (!r) {
342       /* PCReset() does not wipe out the level structure */
343       KSP ksp;
344       PC  pc;
345 
346       ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
347       ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
348       ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr);
349     }
350     /* Cleanup */
351     ierr = VecDestroy(&u);CHKERRQ(ierr);
352     ierr = PetscLogStagePop();CHKERRQ(ierr);
353   }
354   for (r = 1; r <= Nr; ++r) {
355     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
356   }
357   /* Fit convergence rate */
358   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
359   for (f = 0; f < ce->Nf; ++f) {
360     for (r = 0; r <= Nr; ++r) {
361       x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
362       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
363     }
364     ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
365     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
366     alpha[f] = -slope * dim;
367   }
368   ierr = PetscFree2(x, y);CHKERRQ(ierr);
369   ierr = PetscFree2(dm, dof);CHKERRQ(ierr);
370   /* Restore solver */
371   ierr = SNESReset(snes);CHKERRQ(ierr);
372   {
373     /* PCReset() does not wipe out the level structure */
374     KSP ksp;
375     PC  pc;
376 
377     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
378     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
379     ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr);
380     ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */
381   }
382   ierr = SNESSetDM(snes, ce->idm);CHKERRQ(ierr);
383   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
384   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 /*@
389   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
390 
391   Not collective
392 
393   Input Parameter:
394 . ce   - The PetscConvEst object
395 
396   Output Parameter:
397 . alpha - The convergence rate for each field
398 
399   Note: The convergence rate alpha is defined by
400 $ || u_\Delta - u_exact || < C \Delta^alpha
401 where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
402 spatial resolution and \Delta t for the temporal resolution.
403 
404 We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
405 based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
406 
407   Options database keys:
408 + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
409 - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
410 
411   Level: intermediate
412 
413 .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
414 @*/
415 PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
416 {
417   PetscInt       f;
418   PetscErrorCode ierr;
419 
420   PetscFunctionBegin;
421   if (ce->event < 0) {ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);CHKERRQ(ierr);}
422   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
423   ierr = (*ce->ops->getconvrate)(ce, alpha);CHKERRQ(ierr);
424   PetscFunctionReturn(0);
425 }
426 
427 /*@
428   PetscConvEstRateView - Displays the convergence rate to a viewer
429 
430    Collective on SNES
431 
432    Parameter:
433 +  snes - iterative context obtained from SNESCreate()
434 .  alpha - the convergence rate for each field
435 -  viewer - the viewer to display the reason
436 
437    Options Database Keys:
438 .  -snes_convergence_estimate - print the convergence rate
439 
440    Level: developer
441 
442 .seealso: PetscConvEstGetRate()
443 @*/
444 PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
445 {
446   PetscBool      isAscii;
447   PetscErrorCode ierr;
448 
449   PetscFunctionBegin;
450   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
451   if (isAscii) {
452     PetscInt Nf = ce->Nf, f;
453 
454     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
455     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
456     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
457     for (f = 0; f < Nf; ++f) {
458       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
459       ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr);
460     }
461     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
462     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
463     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
464   }
465   PetscFunctionReturn(0);
466 }
467 
468 /*@
469   PetscConvEstCreate - Create a PetscConvEst object
470 
471   Collective
472 
473   Input Parameter:
474 . comm - The communicator for the PetscConvEst object
475 
476   Output Parameter:
477 . ce   - The PetscConvEst object
478 
479   Level: beginner
480 
481 .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
482 @*/
483 PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
484 {
485   PetscErrorCode ierr;
486 
487   PetscFunctionBegin;
488   PetscValidPointer(ce, 2);
489   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
490   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
491   (*ce)->monitor = PETSC_FALSE;
492   (*ce)->r       = 2.0;
493   (*ce)->Nr      = 4;
494   (*ce)->event   = -1;
495   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
496   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
497   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
498   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
499   PetscFunctionReturn(0);
500 }
501