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