xref: /petsc/src/snes/utils/convest.c (revision ca3d3a14de63d3f97e9c05ca58d776b1c45e4686)
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 /*@
166af0ca60SMatthew G. Knepley   PetscConvEstCreate - Create a PetscConvEst object
176af0ca60SMatthew G. Knepley 
186af0ca60SMatthew G. Knepley   Collective on MPI_Comm
196af0ca60SMatthew G. Knepley 
206af0ca60SMatthew G. Knepley   Input Parameter:
216af0ca60SMatthew G. Knepley . comm - The communicator for the PetscConvEst object
226af0ca60SMatthew G. Knepley 
236af0ca60SMatthew G. Knepley   Output Parameter:
246af0ca60SMatthew G. Knepley . ce   - The PetscConvEst object
256af0ca60SMatthew G. Knepley 
266af0ca60SMatthew G. Knepley   Level: beginner
276af0ca60SMatthew G. Knepley 
286af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, create
296af0ca60SMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
306af0ca60SMatthew G. Knepley @*/
316af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
326af0ca60SMatthew G. Knepley {
336af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
346af0ca60SMatthew G. Knepley 
356af0ca60SMatthew G. Knepley   PetscFunctionBegin;
366af0ca60SMatthew G. Knepley   PetscValidPointer(ce, 2);
376af0ca60SMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
386af0ca60SMatthew G. Knepley   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
396af0ca60SMatthew G. Knepley   (*ce)->monitor = PETSC_FALSE;
406af0ca60SMatthew G. Knepley   (*ce)->Nr      = 4;
416af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
426af0ca60SMatthew G. Knepley }
436af0ca60SMatthew G. Knepley 
446af0ca60SMatthew G. Knepley /*@
456af0ca60SMatthew G. Knepley   PetscConvEstDestroy - Destroys a PetscConvEst object
466af0ca60SMatthew G. Knepley 
476af0ca60SMatthew G. Knepley   Collective on PetscConvEst
486af0ca60SMatthew G. Knepley 
496af0ca60SMatthew G. Knepley   Input Parameter:
506af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
516af0ca60SMatthew G. Knepley 
526af0ca60SMatthew G. Knepley   Level: beginner
536af0ca60SMatthew G. Knepley 
546af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, destroy
556af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
566af0ca60SMatthew G. Knepley @*/
576af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
586af0ca60SMatthew G. Knepley {
596af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
606af0ca60SMatthew G. Knepley 
616af0ca60SMatthew G. Knepley   PetscFunctionBegin;
626af0ca60SMatthew G. Knepley   if (!*ce) PetscFunctionReturn(0);
636af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1);
646af0ca60SMatthew G. Knepley   if (--((PetscObject)(*ce))->refct > 0) {
656af0ca60SMatthew G. Knepley     *ce = NULL;
666af0ca60SMatthew G. Knepley     PetscFunctionReturn(0);
676af0ca60SMatthew G. Knepley   }
6895cbbfd3SMatthew G. Knepley   ierr = PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);CHKERRQ(ierr);
696af0ca60SMatthew G. Knepley   ierr = PetscFree((*ce)->errors);CHKERRQ(ierr);
706af0ca60SMatthew G. Knepley   ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr);
716af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
726af0ca60SMatthew G. Knepley }
736af0ca60SMatthew G. Knepley 
746af0ca60SMatthew G. Knepley /*@
756af0ca60SMatthew G. Knepley   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
766af0ca60SMatthew G. Knepley 
776af0ca60SMatthew G. Knepley   Collective on PetscConvEst
786af0ca60SMatthew G. Knepley 
796af0ca60SMatthew G. Knepley   Input Parameters:
806af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
816af0ca60SMatthew G. Knepley 
826af0ca60SMatthew G. Knepley   Level: beginner
836af0ca60SMatthew G. Knepley 
846af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, options
856af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
866af0ca60SMatthew G. Knepley @*/
876af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
886af0ca60SMatthew G. Knepley {
896af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
906af0ca60SMatthew G. Knepley 
916af0ca60SMatthew G. Knepley   PetscFunctionBegin;
926af0ca60SMatthew G. Knepley   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr);
9346079b62SMatthew G. Knepley   ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr);
9446079b62SMatthew G. Knepley   ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr);
956af0ca60SMatthew G. Knepley   ierr = PetscOptionsEnd();
966af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
976af0ca60SMatthew G. Knepley }
986af0ca60SMatthew G. Knepley 
996af0ca60SMatthew G. Knepley /*@
1006af0ca60SMatthew G. Knepley   PetscConvEstView - Views a PetscConvEst object
1016af0ca60SMatthew G. Knepley 
1026af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1036af0ca60SMatthew G. Knepley 
1046af0ca60SMatthew G. Knepley   Input Parameters:
1056af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
1066af0ca60SMatthew G. Knepley - viewer - The PetscViewer object
1076af0ca60SMatthew G. Knepley 
1086af0ca60SMatthew G. Knepley   Level: beginner
1096af0ca60SMatthew G. Knepley 
1106af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, view
1116af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
1126af0ca60SMatthew G. Knepley @*/
1136af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
1146af0ca60SMatthew G. Knepley {
1156af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1166af0ca60SMatthew G. Knepley 
1176af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1186af0ca60SMatthew G. Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr);
1196af0ca60SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr);
1206af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1216af0ca60SMatthew G. Knepley }
1226af0ca60SMatthew G. Knepley 
1236af0ca60SMatthew G. Knepley /*@
1246af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
1256af0ca60SMatthew G. Knepley 
1266af0ca60SMatthew G. Knepley   Not collective
1276af0ca60SMatthew G. Knepley 
1286af0ca60SMatthew G. Knepley   Input Parameter:
1296af0ca60SMatthew G. Knepley . ce   - The PetscConvEst object
1306af0ca60SMatthew G. Knepley 
1316af0ca60SMatthew G. Knepley   Output Parameter:
1326af0ca60SMatthew G. Knepley . snes - The solver
1336af0ca60SMatthew G. Knepley 
1346af0ca60SMatthew G. Knepley   Level: intermediate
1356af0ca60SMatthew G. Knepley 
1366af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence
1376af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1386af0ca60SMatthew G. Knepley @*/
1396af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, SNES *snes)
1406af0ca60SMatthew G. Knepley {
1416af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1426af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
1436af0ca60SMatthew G. Knepley   PetscValidPointer(snes, 2);
1446af0ca60SMatthew G. Knepley   *snes = ce->snes;
1456af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1466af0ca60SMatthew G. Knepley }
1476af0ca60SMatthew G. Knepley 
1486af0ca60SMatthew G. Knepley /*@
1496af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1506af0ca60SMatthew G. Knepley 
1516af0ca60SMatthew G. Knepley   Not collective
1526af0ca60SMatthew G. Knepley 
1536af0ca60SMatthew G. Knepley   Input Parameters:
1546af0ca60SMatthew G. Knepley + ce   - The PetscConvEst object
1556af0ca60SMatthew G. Knepley - snes - The solver
1566af0ca60SMatthew G. Knepley 
1576af0ca60SMatthew G. Knepley   Level: intermediate
1586af0ca60SMatthew G. Knepley 
1596af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1606af0ca60SMatthew G. Knepley 
1616af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence
1626af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1636af0ca60SMatthew G. Knepley @*/
1646af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, SNES snes)
1656af0ca60SMatthew G. Knepley {
1666af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1676af0ca60SMatthew G. Knepley 
1686af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1696af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
1706af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
1716af0ca60SMatthew G. Knepley   ce->snes = snes;
1726af0ca60SMatthew G. Knepley   ierr = SNESGetDM(ce->snes, &ce->idm);CHKERRQ(ierr);
1736af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1746af0ca60SMatthew G. Knepley }
1756af0ca60SMatthew G. Knepley 
1766af0ca60SMatthew G. Knepley /*@
1776af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1786af0ca60SMatthew G. Knepley 
1796af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1806af0ca60SMatthew G. Knepley 
1816af0ca60SMatthew G. Knepley   Input Parameters:
1826af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1836af0ca60SMatthew G. Knepley 
1846af0ca60SMatthew G. Knepley   Level: beginner
1856af0ca60SMatthew G. Knepley 
1866af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, setup
1876af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
1886af0ca60SMatthew G. Knepley @*/
1890955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
1906af0ca60SMatthew G. Knepley {
1916af0ca60SMatthew G. Knepley   PetscDS        prob;
1926af0ca60SMatthew G. Knepley   PetscInt       f;
1936af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1946af0ca60SMatthew G. Knepley 
1956af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1966af0ca60SMatthew G. Knepley   ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr);
1976af0ca60SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &ce->Nf);CHKERRQ(ierr);
1986af0ca60SMatthew G. Knepley   ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
19995cbbfd3SMatthew G. Knepley   ierr = PetscMalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr);
2006af0ca60SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) ce->initGuess[f] = zero_private;
2016af0ca60SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
20295cbbfd3SMatthew G. Knepley     ierr = PetscDSGetExactSolution(prob, f, &ce->exactSol[f], &ce->ctxs[f]);CHKERRQ(ierr);
2036af0ca60SMatthew 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);
2046af0ca60SMatthew G. Knepley   }
2056af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
2066af0ca60SMatthew G. Knepley }
2076af0ca60SMatthew G. Knepley 
2086af0ca60SMatthew G. Knepley /*@
2096af0ca60SMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
2106af0ca60SMatthew G. Knepley 
2116af0ca60SMatthew G. Knepley   Not collective
2126af0ca60SMatthew G. Knepley 
2136af0ca60SMatthew G. Knepley   Input Parameter:
2146af0ca60SMatthew G. Knepley . ce   - The PetscConvEst object
2156af0ca60SMatthew G. Knepley 
2166af0ca60SMatthew G. Knepley   Output Parameter:
21746079b62SMatthew G. Knepley . alpha - The convergence rate for each field
2186af0ca60SMatthew G. Knepley 
2196af0ca60SMatthew G. Knepley   Note: The convergence rate alpha is defined by
2206af0ca60SMatthew G. Knepley $ || u_h - u_exact || < C h^alpha
2216af0ca60SMatthew G. Knepley where u_h is the discrete solution, and h is a measure of the discretization size.
2226af0ca60SMatthew G. Knepley 
2236af0ca60SMatthew G. Knepley We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS,
2246af0ca60SMatthew G. Knepley and then fit the result to our model above using linear regression.
2256af0ca60SMatthew G. Knepley 
2266af0ca60SMatthew G. Knepley   Options database keys:
2276af0ca60SMatthew G. Knepley . -snes_convergence_estimate : Execute convergence estimation and print out the rate
2286af0ca60SMatthew G. Knepley 
2296af0ca60SMatthew G. Knepley   Level: intermediate
2306af0ca60SMatthew G. Knepley 
2316af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence
2326af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
2336af0ca60SMatthew G. Knepley @*/
23446079b62SMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
2356af0ca60SMatthew G. Knepley {
2366af0ca60SMatthew G. Knepley   DM            *dm;
2376af0ca60SMatthew G. Knepley   PetscObject    disc;
2386af0ca60SMatthew G. Knepley   MPI_Comm       comm;
2396af0ca60SMatthew G. Knepley   const char    *uname, *dmname;
2406af0ca60SMatthew G. Knepley   void          *ctx;
2416af0ca60SMatthew G. Knepley   Vec            u;
2426af0ca60SMatthew G. Knepley   PetscReal      t = 0.0, *x, *y, slope, intercept;
24373269098SMatthew G. Knepley   PetscInt      *dof, dim, Nr = ce->Nr, r, f, oldlevel, oldnlev;
244adfa7136SMatthew G. Knepley   PetscLogEvent  event;
2456af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
2466af0ca60SMatthew G. Knepley 
2476af0ca60SMatthew G. Knepley   PetscFunctionBegin;
2486af0ca60SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
2496af0ca60SMatthew G. Knepley   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
2506af0ca60SMatthew G. Knepley   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
2516af0ca60SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
25273269098SMatthew G. Knepley   ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr);
253adfa7136SMatthew G. Knepley   ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr);
2546af0ca60SMatthew G. Knepley   dm[0]  = ce->idm;
25546079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
2566af0ca60SMatthew G. Knepley   /* Loop over meshes */
257adfa7136SMatthew G. Knepley   ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
2586af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
259adfa7136SMatthew G. Knepley     PetscLogStage stage;
260adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
261adfa7136SMatthew G. Knepley 
262adfa7136SMatthew G. Knepley     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
263adfa7136SMatthew G. Knepley     ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);
264adfa7136SMatthew G. Knepley     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
2656af0ca60SMatthew G. Knepley     if (r > 0) {
2666af0ca60SMatthew G. Knepley       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
2676af0ca60SMatthew G. Knepley       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
268e5e52638SMatthew G. Knepley       ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr);
269*ca3d3a14SMatthew G. Knepley       ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr);
2706af0ca60SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
2716af0ca60SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
2728f892730SMatthew G. Knepley       for (f = 0; f <= ce->Nf; ++f) {
2738f892730SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
2748f892730SMatthew G. Knepley         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
2758f892730SMatthew G. Knepley         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
2768f892730SMatthew G. Knepley       }
2776af0ca60SMatthew G. Knepley     }
2786af0ca60SMatthew G. Knepley     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
2796af0ca60SMatthew G. Knepley     /* Create solution */
2806af0ca60SMatthew G. Knepley     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
281e5e52638SMatthew G. Knepley     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
282239a0881SMatthew G. Knepley     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
2836af0ca60SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
2846af0ca60SMatthew G. Knepley     /* Setup solver */
2856af0ca60SMatthew G. Knepley     ierr = SNESReset(ce->snes);CHKERRQ(ierr);
2866af0ca60SMatthew G. Knepley     ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr);
2876af0ca60SMatthew G. Knepley     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
2886af0ca60SMatthew G. Knepley     ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr);
2896af0ca60SMatthew G. Knepley     /* Create initial guess */
29095cbbfd3SMatthew G. Knepley     ierr = DMProjectFunction(dm[r], t, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr);
2916af0ca60SMatthew G. Knepley     ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr);
2925d68e14cSMatthew G. Knepley     ierr = PetscLogEventBegin(event, ce, 0, 0, 0);CHKERRQ(ierr);
29395cbbfd3SMatthew G. Knepley     ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, ce->ctxs, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
294adfa7136SMatthew G. Knepley     ierr = PetscLogEventEnd(event, ce, 0, 0, 0);CHKERRQ(ierr);
295adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
296a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
297a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
298a56d3bf6SMatthew G. Knepley 
299a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
300a56d3bf6SMatthew G. Knepley       ierr = DMGetSection(dm[r], &s);CHKERRQ(ierr);
301a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
302a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr);
303a56d3bf6SMatthew G. Knepley       ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ce->snes));CHKERRQ(ierr);
304a56d3bf6SMatthew G. Knepley       ierr = PetscLogEventSetDof(event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr);
305a56d3bf6SMatthew G. Knepley       ierr = PetscLogEventSetError(event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr);
306adfa7136SMatthew G. Knepley     }
3076af0ca60SMatthew G. Knepley     /* Monitor */
3086af0ca60SMatthew G. Knepley     if (ce->monitor) {
3096af0ca60SMatthew G. Knepley       PetscReal *errors = &ce->errors[r*ce->Nf];
3106af0ca60SMatthew G. Knepley 
311adfa7136SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
312adfa7136SMatthew G. Knepley       if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
3136af0ca60SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3146af0ca60SMatthew G. Knepley         if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
3156af0ca60SMatthew G. Knepley         if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
3168dfade3bSBarry Smith         else                     {ierr = PetscPrintf(comm, "%g", (double)errors[f]);CHKERRQ(ierr);}
3176af0ca60SMatthew G. Knepley       }
318adfa7136SMatthew G. Knepley       if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
319adfa7136SMatthew G. Knepley       ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
3206af0ca60SMatthew G. Knepley     }
32173269098SMatthew G. Knepley     if (!r) {
32273269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
32373269098SMatthew G. Knepley       KSP ksp;
32473269098SMatthew G. Knepley       PC  pc;
32573269098SMatthew G. Knepley 
32673269098SMatthew G. Knepley       ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr);
32773269098SMatthew G. Knepley       ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
32873269098SMatthew G. Knepley       ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr);
32973269098SMatthew G. Knepley     }
3306af0ca60SMatthew G. Knepley     /* Cleanup */
3316af0ca60SMatthew G. Knepley     ierr = VecDestroy(&u);CHKERRQ(ierr);
332adfa7136SMatthew G. Knepley     ierr = PetscLogStagePop();CHKERRQ(ierr);
3336af0ca60SMatthew G. Knepley   }
3346af0ca60SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
3356af0ca60SMatthew G. Knepley     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
3366af0ca60SMatthew G. Knepley   }
3376af0ca60SMatthew G. Knepley   /* Fit convergence rate */
3386af0ca60SMatthew G. Knepley   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
33946079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
3406af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
341adfa7136SMatthew G. Knepley       x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
34246079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
3436af0ca60SMatthew G. Knepley     }
344bebf13c0SMatthew G. Knepley     ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
3456af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
34646079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
34746079b62SMatthew G. Knepley   }
34846079b62SMatthew G. Knepley   ierr = PetscFree2(x, y);CHKERRQ(ierr);
3496af0ca60SMatthew G. Knepley   ierr = PetscFree2(dm, dof);CHKERRQ(ierr);
3502cae373cSMatthew G. Knepley   /* Restore solver */
3512cae373cSMatthew G. Knepley   ierr = SNESReset(ce->snes);CHKERRQ(ierr);
35273269098SMatthew G. Knepley   {
35373269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
35473269098SMatthew G. Knepley     KSP ksp;
35573269098SMatthew G. Knepley     PC  pc;
35673269098SMatthew G. Knepley 
35773269098SMatthew G. Knepley     ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr);
35873269098SMatthew G. Knepley     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
35973269098SMatthew G. Knepley     ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr);
36073269098SMatthew G. Knepley     ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */
36173269098SMatthew G. Knepley   }
3622cae373cSMatthew G. Knepley   ierr = SNESSetDM(ce->snes, ce->idm);CHKERRQ(ierr);
3632cae373cSMatthew G. Knepley   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
3642cae373cSMatthew G. Knepley   ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr);
3656af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
3666af0ca60SMatthew G. Knepley }
3676af0ca60SMatthew G. Knepley 
3686af0ca60SMatthew G. Knepley /*@
3696af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
3706af0ca60SMatthew G. Knepley 
3716af0ca60SMatthew G. Knepley    Collective on SNES
3726af0ca60SMatthew G. Knepley 
3736af0ca60SMatthew G. Knepley    Parameter:
3746af0ca60SMatthew G. Knepley +  snes - iterative context obtained from SNESCreate()
37546079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
3766af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
3776af0ca60SMatthew G. Knepley 
3786af0ca60SMatthew G. Knepley    Options Database Keys:
3796af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
3806af0ca60SMatthew G. Knepley 
3816af0ca60SMatthew G. Knepley    Level: developer
3826af0ca60SMatthew G. Knepley 
3836af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate()
3846af0ca60SMatthew G. Knepley @*/
385a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
3866af0ca60SMatthew G. Knepley {
3876af0ca60SMatthew G. Knepley   PetscBool      isAscii;
3886af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
3896af0ca60SMatthew G. Knepley 
3906af0ca60SMatthew G. Knepley   PetscFunctionBegin;
3916af0ca60SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
3926af0ca60SMatthew G. Knepley   if (isAscii) {
393a56d3bf6SMatthew G. Knepley     DM       dm;
394a56d3bf6SMatthew G. Knepley     PetscInt Nf, f;
395a56d3bf6SMatthew G. Knepley 
396a56d3bf6SMatthew G. Knepley     ierr = SNESGetDM(ce->snes, &dm);CHKERRQ(ierr);
397a56d3bf6SMatthew G. Knepley     ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
3986af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
399adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
400a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
401a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
40246079b62SMatthew G. Knepley       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
4035f9b3039SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr);
40446079b62SMatthew G. Knepley     }
405a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
406adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
4076af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
4086af0ca60SMatthew G. Knepley   }
4096af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4106af0ca60SMatthew G. Knepley }
411