xref: /petsc/src/snes/utils/convest.c (revision 900f6b5bc6e0ed97efa425f9f5f775b415cb9d75)
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);
386af0ca60SMatthew G. Knepley   ierr = PetscFree((*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);
6246079b62SMatthew G. Knepley   ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr);
636af0ca60SMatthew G. Knepley   ierr = PetscOptionsEnd();
646af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
656af0ca60SMatthew G. Knepley }
666af0ca60SMatthew G. Knepley 
676af0ca60SMatthew G. Knepley /*@
686af0ca60SMatthew G. Knepley   PetscConvEstView - Views a PetscConvEst object
696af0ca60SMatthew G. Knepley 
706af0ca60SMatthew G. Knepley   Collective on PetscConvEst
716af0ca60SMatthew G. Knepley 
726af0ca60SMatthew G. Knepley   Input Parameters:
736af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
746af0ca60SMatthew G. Knepley - viewer - The PetscViewer object
756af0ca60SMatthew G. Knepley 
766af0ca60SMatthew G. Knepley   Level: beginner
776af0ca60SMatthew G. Knepley 
786af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
796af0ca60SMatthew G. Knepley @*/
806af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
816af0ca60SMatthew G. Knepley {
826af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
836af0ca60SMatthew G. Knepley 
846af0ca60SMatthew G. Knepley   PetscFunctionBegin;
856af0ca60SMatthew G. Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr);
866af0ca60SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr);
876af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
886af0ca60SMatthew G. Knepley }
896af0ca60SMatthew G. Knepley 
906af0ca60SMatthew G. Knepley /*@
916af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
926af0ca60SMatthew G. Knepley 
936af0ca60SMatthew G. Knepley   Not collective
946af0ca60SMatthew G. Knepley 
956af0ca60SMatthew G. Knepley   Input Parameter:
966af0ca60SMatthew G. Knepley . ce     - The PetscConvEst object
976af0ca60SMatthew G. Knepley 
986af0ca60SMatthew G. Knepley   Output Parameter:
99*900f6b5bSMatthew G. Knepley . solver - The solver
1006af0ca60SMatthew G. Knepley 
1016af0ca60SMatthew G. Knepley   Level: intermediate
1026af0ca60SMatthew G. Knepley 
1036af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1046af0ca60SMatthew G. Knepley @*/
105*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
1066af0ca60SMatthew G. Knepley {
1076af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1086af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
109*900f6b5bSMatthew G. Knepley   PetscValidPointer(solver, 2);
110*900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1116af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1126af0ca60SMatthew G. Knepley }
1136af0ca60SMatthew G. Knepley 
1146af0ca60SMatthew G. Knepley /*@
1156af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1166af0ca60SMatthew G. Knepley 
1176af0ca60SMatthew G. Knepley   Not collective
1186af0ca60SMatthew G. Knepley 
1196af0ca60SMatthew G. Knepley   Input Parameters:
1206af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
121*900f6b5bSMatthew G. Knepley - solver - The solver
1226af0ca60SMatthew G. Knepley 
1236af0ca60SMatthew G. Knepley   Level: intermediate
1246af0ca60SMatthew G. Knepley 
1256af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1266af0ca60SMatthew G. Knepley 
127*900f6b5bSMatthew G. Knepley .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1286af0ca60SMatthew G. Knepley @*/
129*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
1306af0ca60SMatthew G. Knepley {
1316af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1326af0ca60SMatthew G. Knepley 
1336af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1346af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
135*900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
136*900f6b5bSMatthew G. Knepley   ce->solver = solver;
137*900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->setsolver)(ce, solver);CHKERRQ(ierr);
1386af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1396af0ca60SMatthew G. Knepley }
1406af0ca60SMatthew G. Knepley 
1416af0ca60SMatthew G. Knepley /*@
1426af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1436af0ca60SMatthew G. Knepley 
1446af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1456af0ca60SMatthew G. Knepley 
1466af0ca60SMatthew G. Knepley   Input Parameters:
1476af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1486af0ca60SMatthew G. Knepley 
1496af0ca60SMatthew G. Knepley   Level: beginner
1506af0ca60SMatthew G. Knepley 
1516af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
1526af0ca60SMatthew G. Knepley @*/
1530955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
1546af0ca60SMatthew G. Knepley {
155*900f6b5bSMatthew G. Knepley   PetscDS        ds;
156*900f6b5bSMatthew G. Knepley   PetscInt       Nf, f;
1576af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1586af0ca60SMatthew G. Knepley 
1596af0ca60SMatthew G. Knepley   PetscFunctionBegin;
160*900f6b5bSMatthew G. Knepley   ierr = DMGetDS(ce->idm, &ds);CHKERRQ(ierr);
161*900f6b5bSMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
162*900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1636af0ca60SMatthew G. Knepley   ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
16495cbbfd3SMatthew G. Knepley   ierr = PetscMalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr);
165*900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
166*900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
167*900f6b5bSMatthew G. Knepley     ierr = PetscDSGetExactSolution(ds, f, &ce->exactSol[f], &ce->ctxs[f]);CHKERRQ(ierr);
1686af0ca60SMatthew 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);
1696af0ca60SMatthew G. Knepley   }
1706af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1716af0ca60SMatthew G. Knepley }
1726af0ca60SMatthew G. Knepley 
173*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
1746af0ca60SMatthew G. Knepley {
1756af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
1766af0ca60SMatthew G. Knepley 
1776af0ca60SMatthew G. Knepley   PetscFunctionBegin;
178*900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
179*900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
180*900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
181*900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->initguess)(ce, r, dm, u);CHKERRQ(ierr);
182*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
183*900f6b5bSMatthew G. Knepley }
184*900f6b5bSMatthew G. Knepley 
185*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
186*900f6b5bSMatthew G. Knepley {
187*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
188*900f6b5bSMatthew G. Knepley 
189*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
190*900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
191*900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
192*900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
193*900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
194*900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->computeerror)(ce, r, dm, u, errors);CHKERRQ(ierr);
195*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
196*900f6b5bSMatthew G. Knepley }
197*900f6b5bSMatthew G. Knepley 
198*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstMonitor_Private(PetscConvEst ce, PetscInt r)
199*900f6b5bSMatthew G. Knepley {
200*900f6b5bSMatthew G. Knepley   MPI_Comm       comm;
201*900f6b5bSMatthew G. Knepley   PetscInt       f;
202*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
203*900f6b5bSMatthew G. Knepley 
204*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
205*900f6b5bSMatthew G. Knepley   if (ce->monitor) {
206*900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r*ce->Nf];
207*900f6b5bSMatthew G. Knepley 
2086af0ca60SMatthew G. Knepley     ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
209*900f6b5bSMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
210*900f6b5bSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
211*900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
212*900f6b5bSMatthew G. Knepley       if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
213*900f6b5bSMatthew G. Knepley       if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
214*900f6b5bSMatthew G. Knepley       else                     {ierr = PetscPrintf(comm, "%g", (double) errors[f]);CHKERRQ(ierr);}
215*900f6b5bSMatthew G. Knepley     }
216*900f6b5bSMatthew G. Knepley     if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
217*900f6b5bSMatthew G. Knepley     ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
218*900f6b5bSMatthew G. Knepley   }
219*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
220*900f6b5bSMatthew G. Knepley }
221*900f6b5bSMatthew G. Knepley 
222*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
223*900f6b5bSMatthew G. Knepley {
224*900f6b5bSMatthew G. Knepley   PetscClassId   id;
225*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
226*900f6b5bSMatthew G. Knepley 
227*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
228*900f6b5bSMatthew G. Knepley   ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr);
229*900f6b5bSMatthew G. Knepley   if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
230*900f6b5bSMatthew G. Knepley   ierr = SNESGetDM((SNES) ce->solver, &ce->idm);CHKERRQ(ierr);
231*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
232*900f6b5bSMatthew G. Knepley }
233*900f6b5bSMatthew G. Knepley 
234*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
235*900f6b5bSMatthew G. Knepley {
236*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
237*900f6b5bSMatthew G. Knepley 
238*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
239*900f6b5bSMatthew G. Knepley   ierr = DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr);
240*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
241*900f6b5bSMatthew G. Knepley }
242*900f6b5bSMatthew G. Knepley 
243*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
244*900f6b5bSMatthew G. Knepley {
245*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
246*900f6b5bSMatthew G. Knepley 
247*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
248*900f6b5bSMatthew G. Knepley   ierr = DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr);
249*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
250*900f6b5bSMatthew G. Knepley }
251*900f6b5bSMatthew G. Knepley 
252*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
253*900f6b5bSMatthew G. Knepley {
254*900f6b5bSMatthew G. Knepley   SNES           snes = (SNES) ce->solver;
255*900f6b5bSMatthew G. Knepley   DM            *dm;
256*900f6b5bSMatthew G. Knepley   PetscObject    disc;
257*900f6b5bSMatthew G. Knepley   PetscReal     *x, *y, slope, intercept;
258*900f6b5bSMatthew G. Knepley   PetscInt      *dof, Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
259*900f6b5bSMatthew G. Knepley   void          *ctx;
260*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
261*900f6b5bSMatthew G. Knepley 
262*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2636af0ca60SMatthew G. Knepley   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
2646af0ca60SMatthew G. Knepley   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
2656af0ca60SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
26673269098SMatthew G. Knepley   ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr);
267adfa7136SMatthew G. Knepley   ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr);
2686af0ca60SMatthew G. Knepley   /* Loop over meshes */
269*900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
2706af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
271*900f6b5bSMatthew G. Knepley     Vec           u;
272adfa7136SMatthew G. Knepley     PetscLogStage stage;
273adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
274*900f6b5bSMatthew G. Knepley     const char   *dmname, *uname;
275adfa7136SMatthew G. Knepley 
276adfa7136SMatthew G. Knepley     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
277adfa7136SMatthew G. Knepley     ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);
278adfa7136SMatthew G. Knepley     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
2796af0ca60SMatthew G. Knepley     if (r > 0) {
2806af0ca60SMatthew G. Knepley       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
2816af0ca60SMatthew G. Knepley       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
282e5e52638SMatthew G. Knepley       ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr);
283ca3d3a14SMatthew G. Knepley       ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr);
2846af0ca60SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
2856af0ca60SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
2868f892730SMatthew G. Knepley       for (f = 0; f <= ce->Nf; ++f) {
2878f892730SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
288*900f6b5bSMatthew G. Knepley 
2898f892730SMatthew G. Knepley         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
2908f892730SMatthew G. Knepley         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
2918f892730SMatthew G. Knepley       }
2926af0ca60SMatthew G. Knepley     }
2936af0ca60SMatthew G. Knepley     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
2946af0ca60SMatthew G. Knepley     /* Create solution */
2956af0ca60SMatthew G. Knepley     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
296e5e52638SMatthew G. Knepley     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
297239a0881SMatthew G. Knepley     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
2986af0ca60SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
2996af0ca60SMatthew G. Knepley     /* Setup solver */
300*900f6b5bSMatthew G. Knepley     ierr = SNESReset(snes);CHKERRQ(ierr);
301*900f6b5bSMatthew G. Knepley     ierr = SNESSetDM(snes, dm[r]);CHKERRQ(ierr);
3026af0ca60SMatthew G. Knepley     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
303*900f6b5bSMatthew G. Knepley     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
3046af0ca60SMatthew G. Knepley     /* Create initial guess */
305*900f6b5bSMatthew G. Knepley     ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr);
306*900f6b5bSMatthew G. Knepley     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
307*900f6b5bSMatthew G. Knepley     ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
308*900f6b5bSMatthew G. Knepley     ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
309*900f6b5bSMatthew G. Knepley     ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr);
310adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
311a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
312a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
313a56d3bf6SMatthew G. Knepley 
314a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
31592fd8e1eSJed Brown       ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr);
316a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
317a56d3bf6SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr);
318*900f6b5bSMatthew G. Knepley       ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));CHKERRQ(ierr);
319*900f6b5bSMatthew G. Knepley       ierr = PetscLogEventSetDof(ce->event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr);
320*900f6b5bSMatthew G. Knepley       ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr);
321adfa7136SMatthew G. Knepley     }
3226af0ca60SMatthew G. Knepley     /* Monitor */
323*900f6b5bSMatthew G. Knepley     ierr = PetscConvEstMonitor_Private(ce, r);CHKERRQ(ierr);
32473269098SMatthew G. Knepley     if (!r) {
32573269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
32673269098SMatthew G. Knepley       KSP ksp;
32773269098SMatthew G. Knepley       PC  pc;
32873269098SMatthew G. Knepley 
329*900f6b5bSMatthew G. Knepley       ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
33073269098SMatthew G. Knepley       ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
33173269098SMatthew G. Knepley       ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr);
33273269098SMatthew G. Knepley     }
3336af0ca60SMatthew G. Knepley     /* Cleanup */
3346af0ca60SMatthew G. Knepley     ierr = VecDestroy(&u);CHKERRQ(ierr);
335adfa7136SMatthew G. Knepley     ierr = PetscLogStagePop();CHKERRQ(ierr);
3366af0ca60SMatthew G. Knepley   }
3376af0ca60SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
3386af0ca60SMatthew G. Knepley     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
3396af0ca60SMatthew G. Knepley   }
3406af0ca60SMatthew G. Knepley   /* Fit convergence rate */
3416af0ca60SMatthew G. Knepley   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
34246079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
3436af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
344adfa7136SMatthew G. Knepley       x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
34546079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
3466af0ca60SMatthew G. Knepley     }
347bebf13c0SMatthew G. Knepley     ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
3486af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
34946079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
35046079b62SMatthew G. Knepley   }
35146079b62SMatthew G. Knepley   ierr = PetscFree2(x, y);CHKERRQ(ierr);
3526af0ca60SMatthew G. Knepley   ierr = PetscFree2(dm, dof);CHKERRQ(ierr);
3532cae373cSMatthew G. Knepley   /* Restore solver */
354*900f6b5bSMatthew G. Knepley   ierr = SNESReset(snes);CHKERRQ(ierr);
35573269098SMatthew G. Knepley   {
35673269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
35773269098SMatthew G. Knepley     KSP ksp;
35873269098SMatthew G. Knepley     PC  pc;
35973269098SMatthew G. Knepley 
360*900f6b5bSMatthew G. Knepley     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
36173269098SMatthew G. Knepley     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
36273269098SMatthew G. Knepley     ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr);
36373269098SMatthew G. Knepley     ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */
36473269098SMatthew G. Knepley   }
365*900f6b5bSMatthew G. Knepley   ierr = SNESSetDM(snes, ce->idm);CHKERRQ(ierr);
3662cae373cSMatthew G. Knepley   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
367*900f6b5bSMatthew G. Knepley   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
368*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
369*900f6b5bSMatthew G. Knepley }
370*900f6b5bSMatthew G. Knepley 
371*900f6b5bSMatthew G. Knepley /*@
372*900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
373*900f6b5bSMatthew G. Knepley 
374*900f6b5bSMatthew G. Knepley   Not collective
375*900f6b5bSMatthew G. Knepley 
376*900f6b5bSMatthew G. Knepley   Input Parameter:
377*900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
378*900f6b5bSMatthew G. Knepley 
379*900f6b5bSMatthew G. Knepley   Output Parameter:
380*900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
381*900f6b5bSMatthew G. Knepley 
382*900f6b5bSMatthew G. Knepley   Note: The convergence rate alpha is defined by
383*900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
384*900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
385*900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution.
386*900f6b5bSMatthew G. Knepley 
387*900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
388*900f6b5bSMatthew G. Knepley based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
389*900f6b5bSMatthew G. Knepley 
390*900f6b5bSMatthew G. Knepley   Options database keys:
391*900f6b5bSMatthew G. Knepley + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
392*900f6b5bSMatthew G. Knepley - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
393*900f6b5bSMatthew G. Knepley 
394*900f6b5bSMatthew G. Knepley   Level: intermediate
395*900f6b5bSMatthew G. Knepley 
396*900f6b5bSMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
397*900f6b5bSMatthew G. Knepley @*/
398*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
399*900f6b5bSMatthew G. Knepley {
400*900f6b5bSMatthew G. Knepley   PetscInt       f;
401*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
402*900f6b5bSMatthew G. Knepley 
403*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
404*900f6b5bSMatthew G. Knepley   if (ce->event < 0) {ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);CHKERRQ(ierr);}
405*900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
406*900f6b5bSMatthew G. Knepley   ierr = (*ce->ops->getconvrate)(ce, alpha);CHKERRQ(ierr);
4076af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4086af0ca60SMatthew G. Knepley }
4096af0ca60SMatthew G. Knepley 
4106af0ca60SMatthew G. Knepley /*@
4116af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4126af0ca60SMatthew G. Knepley 
4136af0ca60SMatthew G. Knepley    Collective on SNES
4146af0ca60SMatthew G. Knepley 
4156af0ca60SMatthew G. Knepley    Parameter:
4166af0ca60SMatthew G. Knepley +  snes - iterative context obtained from SNESCreate()
41746079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4186af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4196af0ca60SMatthew G. Knepley 
4206af0ca60SMatthew G. Knepley    Options Database Keys:
4216af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4226af0ca60SMatthew G. Knepley 
4236af0ca60SMatthew G. Knepley    Level: developer
4246af0ca60SMatthew G. Knepley 
4256af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate()
4266af0ca60SMatthew G. Knepley @*/
427a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
4286af0ca60SMatthew G. Knepley {
4296af0ca60SMatthew G. Knepley   PetscBool      isAscii;
4306af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
4316af0ca60SMatthew G. Knepley 
4326af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4336af0ca60SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
4346af0ca60SMatthew G. Knepley   if (isAscii) {
435*900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
436a56d3bf6SMatthew G. Knepley 
4376af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
438adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
439a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
440a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
44146079b62SMatthew G. Knepley       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
4425f9b3039SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr);
44346079b62SMatthew G. Knepley     }
444a56d3bf6SMatthew G. Knepley     if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
445adfa7136SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
4466af0ca60SMatthew G. Knepley     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
4476af0ca60SMatthew G. Knepley   }
4486af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4496af0ca60SMatthew G. Knepley }
450*900f6b5bSMatthew G. Knepley 
451*900f6b5bSMatthew G. Knepley /*@
452*900f6b5bSMatthew G. Knepley   PetscConvEstCreate - Create a PetscConvEst object
453*900f6b5bSMatthew G. Knepley 
454*900f6b5bSMatthew G. Knepley   Collective
455*900f6b5bSMatthew G. Knepley 
456*900f6b5bSMatthew G. Knepley   Input Parameter:
457*900f6b5bSMatthew G. Knepley . comm - The communicator for the PetscConvEst object
458*900f6b5bSMatthew G. Knepley 
459*900f6b5bSMatthew G. Knepley   Output Parameter:
460*900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
461*900f6b5bSMatthew G. Knepley 
462*900f6b5bSMatthew G. Knepley   Level: beginner
463*900f6b5bSMatthew G. Knepley 
464*900f6b5bSMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
465*900f6b5bSMatthew G. Knepley @*/
466*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
467*900f6b5bSMatthew G. Knepley {
468*900f6b5bSMatthew G. Knepley   PetscErrorCode ierr;
469*900f6b5bSMatthew G. Knepley 
470*900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
471*900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
472*900f6b5bSMatthew G. Knepley   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
473*900f6b5bSMatthew G. Knepley   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
474*900f6b5bSMatthew G. Knepley   (*ce)->monitor = PETSC_FALSE;
475*900f6b5bSMatthew G. Knepley   (*ce)->Nr      = 4;
476*900f6b5bSMatthew G. Knepley   (*ce)->event   = -1;
477*900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
478*900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
479*900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
480*900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
481*900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
482*900f6b5bSMatthew G. Knepley }
483