xref: /petsc/src/snes/utils/convest.c (revision 608e5a7a8a35b40c4e75a0aae313cbdeb3df9fee)
16af0ca60SMatthew G. Knepley #include <petscconvest.h>            /*I "petscconvest.h" I*/
26af0ca60SMatthew G. Knepley #include <petscdmplex.h>
36af0ca60SMatthew G. Knepley #include <petscds.h>
46af0ca60SMatthew G. Knepley 
56af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h>
66af0ca60SMatthew G. Knepley 
76af0ca60SMatthew G. Knepley static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
86af0ca60SMatthew G. Knepley {
96af0ca60SMatthew G. Knepley   PetscInt c;
106af0ca60SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
116af0ca60SMatthew G. Knepley   return 0;
126af0ca60SMatthew G. Knepley }
136af0ca60SMatthew G. Knepley 
146af0ca60SMatthew G. Knepley /*@
156af0ca60SMatthew G. Knepley   PetscConvEstDestroy - Destroys a PetscConvEst object
166af0ca60SMatthew G. Knepley 
176af0ca60SMatthew G. Knepley   Collective on PetscConvEst
186af0ca60SMatthew G. Knepley 
196af0ca60SMatthew G. Knepley   Input Parameter:
206af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
216af0ca60SMatthew G. Knepley 
226af0ca60SMatthew G. Knepley   Level: beginner
236af0ca60SMatthew G. Knepley 
246af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
256af0ca60SMatthew G. Knepley @*/
266af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
276af0ca60SMatthew G. Knepley {
286af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
296af0ca60SMatthew G. Knepley 
306af0ca60SMatthew G. Knepley   PetscFunctionBegin;
316af0ca60SMatthew G. Knepley   if (!*ce) PetscFunctionReturn(0);
326af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1);
336af0ca60SMatthew G. Knepley   if (--((PetscObject)(*ce))->refct > 0) {
346af0ca60SMatthew G. Knepley     *ce = NULL;
356af0ca60SMatthew G. Knepley     PetscFunctionReturn(0);
366af0ca60SMatthew G. Knepley   }
3795cbbfd3SMatthew G. Knepley   ierr = PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);CHKERRQ(ierr);
387809adefSMatthew G. Knepley   ierr = PetscFree2((*ce)->dofs, (*ce)->errors);CHKERRQ(ierr);
396af0ca60SMatthew G. Knepley   ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr);
406af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
416af0ca60SMatthew G. Knepley }
426af0ca60SMatthew G. Knepley 
436af0ca60SMatthew G. Knepley /*@
446af0ca60SMatthew G. Knepley   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
456af0ca60SMatthew G. Knepley 
466af0ca60SMatthew G. Knepley   Collective on PetscConvEst
476af0ca60SMatthew G. Knepley 
486af0ca60SMatthew G. Knepley   Input Parameters:
496af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
506af0ca60SMatthew G. Knepley 
516af0ca60SMatthew G. Knepley   Level: beginner
526af0ca60SMatthew G. Knepley 
536af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
546af0ca60SMatthew G. Knepley @*/
556af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
566af0ca60SMatthew G. Knepley {
576af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
586af0ca60SMatthew G. Knepley 
596af0ca60SMatthew G. Knepley   PetscFunctionBegin;
606af0ca60SMatthew G. Knepley   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr);
6146079b62SMatthew G. Knepley   ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr);
622e61be88SMatthew G. Knepley   ierr = PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL);CHKERRQ(ierr);
6346079b62SMatthew G. Knepley   ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr);
646af0ca60SMatthew G. Knepley   ierr = PetscOptionsEnd();
656af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
666af0ca60SMatthew G. Knepley }
676af0ca60SMatthew G. Knepley 
686af0ca60SMatthew G. Knepley /*@
696af0ca60SMatthew G. Knepley   PetscConvEstView - Views a PetscConvEst object
706af0ca60SMatthew G. Knepley 
716af0ca60SMatthew G. Knepley   Collective on PetscConvEst
726af0ca60SMatthew G. Knepley 
736af0ca60SMatthew G. Knepley   Input Parameters:
746af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
756af0ca60SMatthew G. Knepley - viewer - The PetscViewer object
766af0ca60SMatthew G. Knepley 
776af0ca60SMatthew G. Knepley   Level: beginner
786af0ca60SMatthew G. Knepley 
796af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
806af0ca60SMatthew G. Knepley @*/
816af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
826af0ca60SMatthew G. Knepley {
836af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
846af0ca60SMatthew G. Knepley 
856af0ca60SMatthew G. Knepley   PetscFunctionBegin;
866af0ca60SMatthew G. Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr);
876af0ca60SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr);
886af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
896af0ca60SMatthew G. Knepley }
906af0ca60SMatthew G. Knepley 
916af0ca60SMatthew G. Knepley /*@
926af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
936af0ca60SMatthew G. Knepley 
946af0ca60SMatthew G. Knepley   Not collective
956af0ca60SMatthew G. Knepley 
966af0ca60SMatthew G. Knepley   Input Parameter:
976af0ca60SMatthew G. Knepley . ce     - The PetscConvEst object
986af0ca60SMatthew G. Knepley 
996af0ca60SMatthew G. Knepley   Output Parameter:
100900f6b5bSMatthew G. Knepley . solver - The solver
1016af0ca60SMatthew G. Knepley 
1026af0ca60SMatthew G. Knepley   Level: intermediate
1036af0ca60SMatthew G. Knepley 
1046af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1056af0ca60SMatthew G. Knepley @*/
106900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
1076af0ca60SMatthew G. Knepley {
1086af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1096af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
110900f6b5bSMatthew G. Knepley   PetscValidPointer(solver, 2);
111900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1126af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1136af0ca60SMatthew G. Knepley }
1146af0ca60SMatthew G. Knepley 
1156af0ca60SMatthew G. Knepley /*@
1166af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1176af0ca60SMatthew G. Knepley 
1186af0ca60SMatthew G. Knepley   Not collective
1196af0ca60SMatthew G. Knepley 
1206af0ca60SMatthew G. Knepley   Input Parameters:
1216af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
122900f6b5bSMatthew G. Knepley - solver - The solver
1236af0ca60SMatthew G. Knepley 
1246af0ca60SMatthew G. Knepley   Level: intermediate
1256af0ca60SMatthew G. Knepley 
1266af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1276af0ca60SMatthew G. Knepley 
128900f6b5bSMatthew G. Knepley .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1296af0ca60SMatthew G. Knepley @*/
130900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
1316af0ca60SMatthew G. Knepley {
1326af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1336af0ca60SMatthew G. Knepley 
1346af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1356af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
136900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
137900f6b5bSMatthew G. Knepley   ce->solver = solver;
138900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->setsolver)(ce, solver);CHKERRQ(ierr);
1396af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1406af0ca60SMatthew G. Knepley }
1416af0ca60SMatthew G. Knepley 
1426af0ca60SMatthew G. Knepley /*@
1436af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1446af0ca60SMatthew G. Knepley 
1456af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1466af0ca60SMatthew G. Knepley 
1476af0ca60SMatthew G. Knepley   Input Parameters:
1486af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1496af0ca60SMatthew G. Knepley 
1506af0ca60SMatthew G. Knepley   Level: beginner
1516af0ca60SMatthew G. Knepley 
1526af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
1536af0ca60SMatthew G. Knepley @*/
1540955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
1556af0ca60SMatthew G. Knepley {
156083401c6SMatthew G. Knepley   PetscInt       Nf, f, Nds, s;
1576af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1586af0ca60SMatthew G. Knepley 
1596af0ca60SMatthew G. Knepley   PetscFunctionBegin;
160083401c6SMatthew G. Knepley   ierr = DMGetNumFields(ce->idm, &Nf);CHKERRQ(ierr);
161900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1627809adefSMatthew G. Knepley   ierr = PetscMalloc2((ce->Nr+1)*ce->Nf, &ce->dofs, (ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
163083401c6SMatthew G. Knepley   ierr = PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr);
164900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
165083401c6SMatthew G. Knepley   ierr = DMGetNumDS(ce->idm, &Nds);CHKERRQ(ierr);
166083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
167083401c6SMatthew G. Knepley     PetscDS         ds;
168083401c6SMatthew G. Knepley     DMLabel         label;
169083401c6SMatthew G. Knepley     IS              fieldIS;
170083401c6SMatthew G. Knepley     const PetscInt *fields;
171083401c6SMatthew G. Knepley     PetscInt        dsNf;
172083401c6SMatthew G. Knepley 
173083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
174083401c6SMatthew G. Knepley     ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
175083401c6SMatthew G. Knepley     if (fieldIS) {ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);}
176083401c6SMatthew G. Knepley     for (f = 0; f < dsNf; ++f) {
177083401c6SMatthew G. Knepley       const PetscInt field = fields[f];
178083401c6SMatthew G. Knepley       ierr = PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]);CHKERRQ(ierr);
179083401c6SMatthew G. Knepley     }
180083401c6SMatthew G. Knepley     if (fieldIS) {ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);}
181083401c6SMatthew G. Knepley   }
182900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1836af0ca60SMatthew G. Knepley     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);
1846af0ca60SMatthew G. Knepley   }
1856af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1866af0ca60SMatthew G. Knepley }
1876af0ca60SMatthew G. Knepley 
188900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
1896af0ca60SMatthew G. Knepley {
1906af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1916af0ca60SMatthew G. Knepley 
1926af0ca60SMatthew G. Knepley   PetscFunctionBegin;
193900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
194900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
195900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
196900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->initguess)(ce, r, dm, u);CHKERRQ(ierr);
197900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
198900f6b5bSMatthew G. Knepley }
199900f6b5bSMatthew G. Knepley 
200900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
201900f6b5bSMatthew G. Knepley {
202900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
203900f6b5bSMatthew G. Knepley 
204900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
205900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
206900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
207900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
208900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
209900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->computeerror)(ce, r, dm, u, errors);CHKERRQ(ierr);
210900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
211900f6b5bSMatthew G. Knepley }
212900f6b5bSMatthew G. Knepley 
213f2cacb80SMatthew G. Knepley /*@
214f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
215f2cacb80SMatthew G. Knepley 
216f2cacb80SMatthew G. Knepley   Collective on PetscConvEst
217f2cacb80SMatthew G. Knepley 
218f2cacb80SMatthew G. Knepley   Input Parameter:
219f2cacb80SMatthew G. Knepley + ce - The PetscConvEst object
220f2cacb80SMatthew G. Knepley - r  - The refinement level
221f2cacb80SMatthew G. Knepley 
222f2cacb80SMatthew G. Knepley   Options database keys:
223ee300463SSatish Balay . -convest_monitor - Activate the monitor
224f2cacb80SMatthew G. Knepley 
225f2cacb80SMatthew G. Knepley   Level: intermediate
226f2cacb80SMatthew G. Knepley 
227f2cacb80SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
228f2cacb80SMatthew G. Knepley @*/
229f2cacb80SMatthew G. Knepley PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
230900f6b5bSMatthew G. Knepley {
231900f6b5bSMatthew G. Knepley   MPI_Comm       comm;
232900f6b5bSMatthew G. Knepley   PetscInt       f;
233900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
234900f6b5bSMatthew G. Knepley 
235900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
236900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2377809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r*ce->Nf];
238900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r*ce->Nf];
239900f6b5bSMatthew G. Knepley 
2406af0ca60SMatthew G. Knepley     ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
2417809adefSMatthew G. Knepley     ierr = PetscPrintf(comm, "N: ");CHKERRQ(ierr);
2427809adefSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
2437809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2447809adefSMatthew G. Knepley       if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
2457809adefSMatthew G. Knepley       ierr = PetscPrintf(comm, "%7D", dofs[f]);CHKERRQ(ierr);
2467809adefSMatthew G. Knepley     }
2477809adefSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
2487809adefSMatthew G. Knepley     ierr = PetscPrintf(comm, "  ");CHKERRQ(ierr);
249900f6b5bSMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
250900f6b5bSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
251900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
252900f6b5bSMatthew G. Knepley       if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
253900f6b5bSMatthew G. Knepley       if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
254900f6b5bSMatthew G. Knepley       else                     {ierr = PetscPrintf(comm, "%g", (double) errors[f]);CHKERRQ(ierr);}
255900f6b5bSMatthew G. Knepley     }
256900f6b5bSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
257900f6b5bSMatthew G. Knepley     ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
258900f6b5bSMatthew G. Knepley   }
259900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
260900f6b5bSMatthew G. Knepley }
261900f6b5bSMatthew G. Knepley 
262900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
263900f6b5bSMatthew G. Knepley {
264900f6b5bSMatthew G. Knepley   PetscClassId   id;
265900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
266900f6b5bSMatthew G. Knepley 
267900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
268900f6b5bSMatthew G. Knepley   ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr);
269900f6b5bSMatthew G. Knepley   if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
270900f6b5bSMatthew G. Knepley   ierr = SNESGetDM((SNES) ce->solver, &ce->idm);CHKERRQ(ierr);
271900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
272900f6b5bSMatthew G. Knepley }
273900f6b5bSMatthew G. Knepley 
274900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
275900f6b5bSMatthew G. Knepley {
276900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
277900f6b5bSMatthew G. Knepley 
278900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
279900f6b5bSMatthew G. Knepley   ierr = DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr);
280900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
281900f6b5bSMatthew G. Knepley }
282900f6b5bSMatthew G. Knepley 
283900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
284900f6b5bSMatthew G. Knepley {
285900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
286900f6b5bSMatthew G. Knepley 
287900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
288900f6b5bSMatthew G. Knepley   ierr = DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr);
289900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
290900f6b5bSMatthew G. Knepley }
291900f6b5bSMatthew G. Knepley 
292478db826SMatthew G. Knepley static PetscErrorCode PetscConvEstSetJacobianNullspace_Private(PetscConvEst ce, SNES snes)
293478db826SMatthew G. Knepley {
294478db826SMatthew G. Knepley   DM             dm;
295478db826SMatthew G. Knepley   PetscInt       f;
296478db826SMatthew G. Knepley   PetscErrorCode ierr;
297478db826SMatthew G. Knepley 
298478db826SMatthew G. Knepley   PetscFunctionBegin;
299478db826SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
300478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
301478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
302478db826SMatthew G. Knepley 
303478db826SMatthew G. Knepley     ierr = DMGetNullSpaceConstructor(dm, f, &nspconstr);CHKERRQ(ierr);
304478db826SMatthew G. Knepley     if (nspconstr) {
305478db826SMatthew G. Knepley       MatNullSpace nullsp;
306478db826SMatthew G. Knepley       Mat          J;
307478db826SMatthew G. Knepley 
308478db826SMatthew G. Knepley       ierr = (*nspconstr)(dm, f, f,&nullsp);CHKERRQ(ierr);
309478db826SMatthew G. Knepley       ierr = SNESSetUp(snes);CHKERRQ(ierr);
310478db826SMatthew G. Knepley       ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
311478db826SMatthew G. Knepley       ierr = MatSetNullSpace(J, nullsp);CHKERRQ(ierr);
312478db826SMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
313478db826SMatthew G. Knepley       break;
314478db826SMatthew G. Knepley     }
315478db826SMatthew G. Knepley   }
316478db826SMatthew G. Knepley   PetscFunctionReturn(0);
317478db826SMatthew G. Knepley }
318478db826SMatthew G. Knepley 
319900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
320900f6b5bSMatthew G. Knepley {
321900f6b5bSMatthew G. Knepley   SNES           snes = (SNES) ce->solver;
322900f6b5bSMatthew G. Knepley   DM            *dm;
323900f6b5bSMatthew G. Knepley   PetscObject    disc;
324900f6b5bSMatthew G. Knepley   PetscReal     *x, *y, slope, intercept;
3257809adefSMatthew G. Knepley   PetscInt       Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
326900f6b5bSMatthew G. Knepley   void          *ctx;
327900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
328900f6b5bSMatthew G. Knepley 
329900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
3302e61be88SMatthew G. Knepley   if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
3316af0ca60SMatthew G. Knepley   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
3326af0ca60SMatthew G. Knepley   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
3336af0ca60SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
33473269098SMatthew G. Knepley   ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr);
3357809adefSMatthew G. Knepley   ierr = PetscMalloc1((Nr+1), &dm);CHKERRQ(ierr);
3366af0ca60SMatthew G. Knepley   /* Loop over meshes */
337900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3386af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
339900f6b5bSMatthew G. Knepley     Vec           u;
340e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG)
341adfa7136SMatthew G. Knepley     PetscLogStage stage;
342e5ed2c37SJose E. Roman #endif
343adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
344900f6b5bSMatthew G. Knepley     const char   *dmname, *uname;
345adfa7136SMatthew G. Knepley 
346adfa7136SMatthew G. Knepley     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
347*608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
348*608e5a7aSMatthew G. Knepley     ierr = PetscLogStageGetId(stageName, &stage);CHKERRQ(ierr);
349*608e5a7aSMatthew G. Knepley     if (stage < 0) {ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);}
350*608e5a7aSMatthew G. Knepley #endif
351adfa7136SMatthew G. Knepley     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
3526af0ca60SMatthew G. Knepley     if (r > 0) {
3536af0ca60SMatthew G. Knepley       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
3546af0ca60SMatthew G. Knepley       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
355ca3d3a14SMatthew G. Knepley       ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr);
3566af0ca60SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
3576af0ca60SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
358478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3598cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
360900f6b5bSMatthew G. Knepley 
3618f892730SMatthew G. Knepley         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
3628f892730SMatthew G. Knepley         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
3638f892730SMatthew G. Knepley       }
3646af0ca60SMatthew G. Knepley     }
3656af0ca60SMatthew G. Knepley     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
3666af0ca60SMatthew G. Knepley     /* Create solution */
3676af0ca60SMatthew G. Knepley     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
368e5e52638SMatthew G. Knepley     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
369239a0881SMatthew G. Knepley     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
3706af0ca60SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
3716af0ca60SMatthew G. Knepley     /* Setup solver */
372900f6b5bSMatthew G. Knepley     ierr = SNESReset(snes);CHKERRQ(ierr);
373900f6b5bSMatthew G. Knepley     ierr = SNESSetDM(snes, dm[r]);CHKERRQ(ierr);
3746af0ca60SMatthew G. Knepley     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
375900f6b5bSMatthew G. Knepley     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
376478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
377478db826SMatthew G. Knepley     ierr = PetscConvEstSetJacobianNullspace_Private(ce, snes);CHKERRQ(ierr);
3786af0ca60SMatthew G. Knepley     /* Create initial guess */
379900f6b5bSMatthew G. Knepley     ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr);
380900f6b5bSMatthew G. Knepley     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
381900f6b5bSMatthew G. Knepley     ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
382900f6b5bSMatthew G. Knepley     ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
383900f6b5bSMatthew G. Knepley     ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
384adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
385a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
386a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
387a56d3bf6SMatthew G. Knepley 
388a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
38992fd8e1eSJed Brown       ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr);
390a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
391a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr);
392ffc4695bSBarry Smith       ierr = MPI_Allreduce(&lsize, &ce->dofs[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));CHKERRMPI(ierr);
3937809adefSMatthew G. Knepley       ierr = PetscLogEventSetDof(ce->event, f, ce->dofs[r*ce->Nf+f]);CHKERRQ(ierr);
394900f6b5bSMatthew G. Knepley       ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr);
395adfa7136SMatthew G. Knepley     }
3966af0ca60SMatthew G. Knepley     /* Monitor */
397f2cacb80SMatthew G. Knepley     ierr = PetscConvEstMonitorDefault(ce, r);CHKERRQ(ierr);
39873269098SMatthew G. Knepley     if (!r) {
39973269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
40073269098SMatthew G. Knepley       KSP ksp;
40173269098SMatthew G. Knepley       PC  pc;
40273269098SMatthew G. Knepley 
403900f6b5bSMatthew G. Knepley       ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
40473269098SMatthew G. Knepley       ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
40573269098SMatthew G. Knepley       ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr);
40673269098SMatthew G. Knepley     }
4076af0ca60SMatthew G. Knepley     /* Cleanup */
4086af0ca60SMatthew G. Knepley     ierr = VecDestroy(&u);CHKERRQ(ierr);
409adfa7136SMatthew G. Knepley     ierr = PetscLogStagePop();CHKERRQ(ierr);
4106af0ca60SMatthew G. Knepley   }
4116af0ca60SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
4126af0ca60SMatthew G. Knepley     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
4136af0ca60SMatthew G. Knepley   }
4146af0ca60SMatthew G. Knepley   /* Fit convergence rate */
4156af0ca60SMatthew G. Knepley   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
41646079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
4176af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
4187809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r*ce->Nf+f]);
41946079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
4206af0ca60SMatthew G. Knepley     }
421bebf13c0SMatthew G. Knepley     ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
4226af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
42346079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
42446079b62SMatthew G. Knepley   }
42546079b62SMatthew G. Knepley   ierr = PetscFree2(x, y);CHKERRQ(ierr);
4267809adefSMatthew G. Knepley   ierr = PetscFree(dm);CHKERRQ(ierr);
4272cae373cSMatthew G. Knepley   /* Restore solver */
428900f6b5bSMatthew G. Knepley   ierr = SNESReset(snes);CHKERRQ(ierr);
42973269098SMatthew G. Knepley   {
43073269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
43173269098SMatthew G. Knepley     KSP ksp;
43273269098SMatthew G. Knepley     PC  pc;
43373269098SMatthew G. Knepley 
434900f6b5bSMatthew G. Knepley     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
43573269098SMatthew G. Knepley     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
43673269098SMatthew G. Knepley     ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr);
43773269098SMatthew G. Knepley     ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */
43873269098SMatthew G. Knepley   }
439900f6b5bSMatthew G. Knepley   ierr = SNESSetDM(snes, ce->idm);CHKERRQ(ierr);
4402cae373cSMatthew G. Knepley   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
441900f6b5bSMatthew G. Knepley   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
442478db826SMatthew G. Knepley   ierr = PetscConvEstSetJacobianNullspace_Private(ce, snes);CHKERRQ(ierr);
443900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
444900f6b5bSMatthew G. Knepley }
445900f6b5bSMatthew G. Knepley 
446900f6b5bSMatthew G. Knepley /*@
447900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
448900f6b5bSMatthew G. Knepley 
449900f6b5bSMatthew G. Knepley   Not collective
450900f6b5bSMatthew G. Knepley 
451900f6b5bSMatthew G. Knepley   Input Parameter:
452900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
453900f6b5bSMatthew G. Knepley 
454900f6b5bSMatthew G. Knepley   Output Parameter:
455900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
456900f6b5bSMatthew G. Knepley 
457900f6b5bSMatthew G. Knepley   Note: The convergence rate alpha is defined by
458900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
459900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
460900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution.
461900f6b5bSMatthew G. Knepley 
462900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
463900f6b5bSMatthew G. Knepley based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
464900f6b5bSMatthew G. Knepley 
465900f6b5bSMatthew G. Knepley   Options database keys:
466900f6b5bSMatthew G. Knepley + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
467900f6b5bSMatthew G. Knepley - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
468900f6b5bSMatthew G. Knepley 
469900f6b5bSMatthew G. Knepley   Level: intermediate
470900f6b5bSMatthew G. Knepley 
471900f6b5bSMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
472900f6b5bSMatthew G. Knepley @*/
473900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
474900f6b5bSMatthew G. Knepley {
475900f6b5bSMatthew G. Knepley   PetscInt       f;
476900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
477900f6b5bSMatthew G. Knepley 
478900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
479900f6b5bSMatthew G. Knepley   if (ce->event < 0) {ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);CHKERRQ(ierr);}
480900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
481900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->getconvrate)(ce, alpha);CHKERRQ(ierr);
4826af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4836af0ca60SMatthew G. Knepley }
4846af0ca60SMatthew G. Knepley 
4856af0ca60SMatthew G. Knepley /*@
4866af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4876af0ca60SMatthew G. Knepley 
4886af0ca60SMatthew G. Knepley    Collective on SNES
4896af0ca60SMatthew G. Knepley 
4906af0ca60SMatthew G. Knepley    Parameter:
4916af0ca60SMatthew G. Knepley +  snes - iterative context obtained from SNESCreate()
49246079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4936af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4946af0ca60SMatthew G. Knepley 
4956af0ca60SMatthew G. Knepley    Options Database Keys:
4966af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4976af0ca60SMatthew G. Knepley 
4986af0ca60SMatthew G. Knepley    Level: developer
4996af0ca60SMatthew G. Knepley 
5006af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate()
5016af0ca60SMatthew G. Knepley @*/
502a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
5036af0ca60SMatthew G. Knepley {
5046af0ca60SMatthew G. Knepley   PetscBool      isAscii;
5056af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
5066af0ca60SMatthew G. Knepley 
5076af0ca60SMatthew G. Knepley   PetscFunctionBegin;
5086af0ca60SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
5096af0ca60SMatthew G. Knepley   if (isAscii) {
510900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
511a56d3bf6SMatthew G. Knepley 
5126af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
513adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
514a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
515a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
51646079b62SMatthew G. Knepley       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
5175f9b3039SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr);
51846079b62SMatthew G. Knepley     }
519a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
520adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
5216af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
5226af0ca60SMatthew G. Knepley   }
5236af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
5246af0ca60SMatthew G. Knepley }
525900f6b5bSMatthew G. Knepley 
526900f6b5bSMatthew G. Knepley /*@
527900f6b5bSMatthew G. Knepley   PetscConvEstCreate - Create a PetscConvEst object
528900f6b5bSMatthew G. Knepley 
529900f6b5bSMatthew G. Knepley   Collective
530900f6b5bSMatthew G. Knepley 
531900f6b5bSMatthew G. Knepley   Input Parameter:
532900f6b5bSMatthew G. Knepley . comm - The communicator for the PetscConvEst object
533900f6b5bSMatthew G. Knepley 
534900f6b5bSMatthew G. Knepley   Output Parameter:
535900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
536900f6b5bSMatthew G. Knepley 
537900f6b5bSMatthew G. Knepley   Level: beginner
538900f6b5bSMatthew G. Knepley 
539900f6b5bSMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
540900f6b5bSMatthew G. Knepley @*/
541900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
542900f6b5bSMatthew G. Knepley {
543900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
544900f6b5bSMatthew G. Knepley 
545900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
546900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
547900f6b5bSMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
548900f6b5bSMatthew G. Knepley   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
549900f6b5bSMatthew G. Knepley   (*ce)->monitor = PETSC_FALSE;
5502e61be88SMatthew G. Knepley   (*ce)->r       = 2.0;
551900f6b5bSMatthew G. Knepley   (*ce)->Nr      = 4;
552900f6b5bSMatthew G. Knepley   (*ce)->event   = -1;
553900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
554900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
555900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
556900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
557900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
558900f6b5bSMatthew G. Knepley }
559