xref: /petsc/src/snes/utils/convest.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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   PetscFunctionBegin;
296af0ca60SMatthew G. Knepley   if (!*ce) PetscFunctionReturn(0);
306af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1);
316af0ca60SMatthew G. Knepley   if (--((PetscObject)(*ce))->refct > 0) {
326af0ca60SMatthew G. Knepley     *ce = NULL;
336af0ca60SMatthew G. Knepley     PetscFunctionReturn(0);
346af0ca60SMatthew G. Knepley   }
359566063dSJacob Faibussowitsch   PetscCall(PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs));
369566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*ce)->dofs, (*ce)->errors));
379566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(ce));
386af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
396af0ca60SMatthew G. Knepley }
406af0ca60SMatthew G. Knepley 
416af0ca60SMatthew G. Knepley /*@
426af0ca60SMatthew G. Knepley   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
436af0ca60SMatthew G. Knepley 
446af0ca60SMatthew G. Knepley   Collective on PetscConvEst
456af0ca60SMatthew G. Knepley 
466af0ca60SMatthew G. Knepley   Input Parameters:
476af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
486af0ca60SMatthew G. Knepley 
496af0ca60SMatthew G. Knepley   Level: beginner
506af0ca60SMatthew G. Knepley 
516af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
526af0ca60SMatthew G. Knepley @*/
536af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
546af0ca60SMatthew G. Knepley {
556af0ca60SMatthew G. Knepley   PetscErrorCode ierr;
566af0ca60SMatthew G. Knepley 
576af0ca60SMatthew G. Knepley   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");PetscCall(ierr);
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL));
639566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
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   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer));
849566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1));
856af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
866af0ca60SMatthew G. Knepley }
876af0ca60SMatthew G. Knepley 
886af0ca60SMatthew G. Knepley /*@
896af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
906af0ca60SMatthew G. Knepley 
916af0ca60SMatthew G. Knepley   Not collective
926af0ca60SMatthew G. Knepley 
936af0ca60SMatthew G. Knepley   Input Parameter:
946af0ca60SMatthew G. Knepley . ce     - The PetscConvEst object
956af0ca60SMatthew G. Knepley 
966af0ca60SMatthew G. Knepley   Output Parameter:
97900f6b5bSMatthew G. Knepley . solver - The solver
986af0ca60SMatthew G. Knepley 
996af0ca60SMatthew G. Knepley   Level: intermediate
1006af0ca60SMatthew G. Knepley 
1016af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1026af0ca60SMatthew G. Knepley @*/
103900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
1046af0ca60SMatthew G. Knepley {
1056af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1066af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
107900f6b5bSMatthew G. Knepley   PetscValidPointer(solver, 2);
108900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1096af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1106af0ca60SMatthew G. Knepley }
1116af0ca60SMatthew G. Knepley 
1126af0ca60SMatthew G. Knepley /*@
1136af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1146af0ca60SMatthew G. Knepley 
1156af0ca60SMatthew G. Knepley   Not collective
1166af0ca60SMatthew G. Knepley 
1176af0ca60SMatthew G. Knepley   Input Parameters:
1186af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
119900f6b5bSMatthew G. Knepley - solver - The solver
1206af0ca60SMatthew G. Knepley 
1216af0ca60SMatthew G. Knepley   Level: intermediate
1226af0ca60SMatthew G. Knepley 
1236af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1246af0ca60SMatthew G. Knepley 
125900f6b5bSMatthew G. Knepley .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate()
1266af0ca60SMatthew G. Knepley @*/
127900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
1286af0ca60SMatthew G. Knepley {
1296af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1306af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
131900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
132900f6b5bSMatthew G. Knepley   ce->solver = solver;
1339566063dSJacob Faibussowitsch   PetscCall((*ce->ops->setsolver)(ce, solver));
1346af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1356af0ca60SMatthew G. Knepley }
1366af0ca60SMatthew G. Knepley 
1376af0ca60SMatthew G. Knepley /*@
1386af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1396af0ca60SMatthew G. Knepley 
1406af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1416af0ca60SMatthew G. Knepley 
1426af0ca60SMatthew G. Knepley   Input Parameters:
1436af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1446af0ca60SMatthew G. Knepley 
1456af0ca60SMatthew G. Knepley   Level: beginner
1466af0ca60SMatthew G. Knepley 
1476af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
1486af0ca60SMatthew G. Knepley @*/
1490955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
1506af0ca60SMatthew G. Knepley {
151083401c6SMatthew G. Knepley   PetscInt       Nf, f, Nds, s;
1526af0ca60SMatthew G. Knepley 
1536af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1549566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(ce->idm, &Nf));
155900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2((ce->Nr+1)*ce->Nf, &ce->dofs, (ce->Nr+1)*ce->Nf, &ce->errors));
1579566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs));
158900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
1599566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(ce->idm, &Nds));
160083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
161083401c6SMatthew G. Knepley     PetscDS         ds;
162083401c6SMatthew G. Knepley     DMLabel         label;
163083401c6SMatthew G. Knepley     IS              fieldIS;
164083401c6SMatthew G. Knepley     const PetscInt *fields;
165083401c6SMatthew G. Knepley     PetscInt        dsNf;
166083401c6SMatthew G. Knepley 
1679566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds));
1689566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &dsNf));
1699566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
170083401c6SMatthew G. Knepley     for (f = 0; f < dsNf; ++f) {
171083401c6SMatthew G. Knepley       const PetscInt field = fields[f];
1729566063dSJacob Faibussowitsch       PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]));
173083401c6SMatthew G. Knepley     }
1749566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
175083401c6SMatthew G. Knepley   }
176900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
177*08401ef6SPierre Jolivet     PetscCheck(ce->exactSol[f],PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %D", f);
1786af0ca60SMatthew G. Knepley   }
1796af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1806af0ca60SMatthew G. Knepley }
1816af0ca60SMatthew G. Knepley 
182900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
1836af0ca60SMatthew G. Knepley {
1846af0ca60SMatthew G. Knepley   PetscFunctionBegin;
185900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
186900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
187900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1889566063dSJacob Faibussowitsch   PetscCall((*ce->ops->initguess)(ce, r, dm, u));
189900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
190900f6b5bSMatthew G. Knepley }
191900f6b5bSMatthew G. Knepley 
192900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
193900f6b5bSMatthew G. Knepley {
194900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
195900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
196900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
197900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
198900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
1999566063dSJacob Faibussowitsch   PetscCall((*ce->ops->computeerror)(ce, r, dm, u, errors));
200900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
201900f6b5bSMatthew G. Knepley }
202900f6b5bSMatthew G. Knepley 
203f2cacb80SMatthew G. Knepley /*@
204f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
205f2cacb80SMatthew G. Knepley 
206f2cacb80SMatthew G. Knepley   Collective on PetscConvEst
207f2cacb80SMatthew G. Knepley 
208d8d19677SJose E. Roman   Input Parameters:
209f2cacb80SMatthew G. Knepley + ce - The PetscConvEst object
210f2cacb80SMatthew G. Knepley - r  - The refinement level
211f2cacb80SMatthew G. Knepley 
212f2cacb80SMatthew G. Knepley   Options database keys:
213ee300463SSatish Balay . -convest_monitor - Activate the monitor
214f2cacb80SMatthew G. Knepley 
215f2cacb80SMatthew G. Knepley   Level: intermediate
216f2cacb80SMatthew G. Knepley 
217f2cacb80SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
218f2cacb80SMatthew G. Knepley @*/
219f2cacb80SMatthew G. Knepley PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
220900f6b5bSMatthew G. Knepley {
221900f6b5bSMatthew G. Knepley   MPI_Comm       comm;
222900f6b5bSMatthew G. Knepley   PetscInt       f;
223900f6b5bSMatthew G. Knepley 
224900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
225900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2267809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r*ce->Nf];
227900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r*ce->Nf];
228900f6b5bSMatthew G. Knepley 
2299566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) ce, &comm));
2309566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "N: "));
2319566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
2327809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2339566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2349566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%7D", dofs[f]));
2357809adefSMatthew G. Knepley     }
2369566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2379566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "  "));
2389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Error: "));
2399566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
240900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2419566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2429566063dSJacob Faibussowitsch       if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11"));
2439566063dSJacob Faibussowitsch       else                     PetscCall(PetscPrintf(comm, "%g", (double) errors[f]));
244900f6b5bSMatthew G. Knepley     }
2459566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2469566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n"));
247900f6b5bSMatthew G. Knepley   }
248900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
249900f6b5bSMatthew G. Knepley }
250900f6b5bSMatthew G. Knepley 
251900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
252900f6b5bSMatthew G. Knepley {
253900f6b5bSMatthew G. Knepley   PetscClassId   id;
254900f6b5bSMatthew G. Knepley 
255900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
257*08401ef6SPierre Jolivet   PetscCheck(id == SNES_CLASSID,PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
2589566063dSJacob Faibussowitsch   PetscCall(SNESGetDM((SNES) ce->solver, &ce->idm));
259900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
260900f6b5bSMatthew G. Knepley }
261900f6b5bSMatthew G. Knepley 
262900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
263900f6b5bSMatthew G. Knepley {
264900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u));
266900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
267900f6b5bSMatthew G. Knepley }
268900f6b5bSMatthew G. Knepley 
269900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
270900f6b5bSMatthew G. Knepley {
271900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2729566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors));
273900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
274900f6b5bSMatthew G. Knepley }
275900f6b5bSMatthew G. Knepley 
276478db826SMatthew G. Knepley static PetscErrorCode PetscConvEstSetJacobianNullspace_Private(PetscConvEst ce, SNES snes)
277478db826SMatthew G. Knepley {
278478db826SMatthew G. Knepley   DM             dm;
279478db826SMatthew G. Knepley   PetscInt       f;
280478db826SMatthew G. Knepley 
281478db826SMatthew G. Knepley   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
283478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
284478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
285478db826SMatthew G. Knepley 
2869566063dSJacob Faibussowitsch     PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr));
287478db826SMatthew G. Knepley     if (nspconstr) {
288478db826SMatthew G. Knepley       MatNullSpace nullsp;
289478db826SMatthew G. Knepley       Mat          J;
290478db826SMatthew G. Knepley 
2919566063dSJacob Faibussowitsch       PetscCall((*nspconstr)(dm, f, f,&nullsp));
2929566063dSJacob Faibussowitsch       PetscCall(SNESSetUp(snes));
2939566063dSJacob Faibussowitsch       PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
2949566063dSJacob Faibussowitsch       PetscCall(MatSetNullSpace(J, nullsp));
2959566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullsp));
296478db826SMatthew G. Knepley       break;
297478db826SMatthew G. Knepley     }
298478db826SMatthew G. Knepley   }
299478db826SMatthew G. Knepley   PetscFunctionReturn(0);
300478db826SMatthew G. Knepley }
301478db826SMatthew G. Knepley 
302900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
303900f6b5bSMatthew G. Knepley {
304900f6b5bSMatthew G. Knepley   SNES           snes = (SNES) ce->solver;
305900f6b5bSMatthew G. Knepley   DM            *dm;
306900f6b5bSMatthew G. Knepley   PetscObject    disc;
307900f6b5bSMatthew G. Knepley   PetscReal     *x, *y, slope, intercept;
3087809adefSMatthew G. Knepley   PetscInt       Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
309900f6b5bSMatthew G. Knepley   void          *ctx;
310900f6b5bSMatthew G. Knepley 
311900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
312*08401ef6SPierre Jolivet   PetscCheck(ce->r == 2.0,PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
3139566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
3149566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
3159566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
3169566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
3179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr+1), &dm));
3186af0ca60SMatthew G. Knepley   /* Loop over meshes */
319900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3206af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
321900f6b5bSMatthew G. Knepley     Vec           u;
322e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG)
323adfa7136SMatthew G. Knepley     PetscLogStage stage;
324e5ed2c37SJose E. Roman #endif
325adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
326900f6b5bSMatthew G. Knepley     const char   *dmname, *uname;
327adfa7136SMatthew G. Knepley 
3289566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r));
329608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
3309566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
3319566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
332608e5a7aSMatthew G. Knepley #endif
3339566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3346af0ca60SMatthew G. Knepley     if (r > 0) {
335b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
3369566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]));
3379566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r-1]));
338b2df8587SMatthew G. Knepley       } else {
339b2df8587SMatthew G. Knepley         DM cdm, rcdm;
340b2df8587SMatthew G. Knepley 
3419566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r-1], &dm[r]));
3429566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r-1], dm[r]));
3439566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r-1], &cdm));
3449566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r],   &rcdm));
3459566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
346b2df8587SMatthew G. Knepley       }
3479566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
3489566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) dm[r-1], &dmname));
3499566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) dm[r], dmname));
350478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3518cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
352900f6b5bSMatthew G. Knepley 
3539566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr));
3549566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r],   f,  nspconstr));
3558f892730SMatthew G. Knepley       }
3566af0ca60SMatthew G. Knepley     }
3579566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
3586af0ca60SMatthew G. Knepley     /* Create solution */
3599566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
3609566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
3619566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
3629566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) u, uname));
3636af0ca60SMatthew G. Knepley     /* Setup solver */
3649566063dSJacob Faibussowitsch     PetscCall(SNESReset(snes));
3659566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm[r]));
3669566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx));
3679566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
368478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
3699566063dSJacob Faibussowitsch     PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
3706af0ca60SMatthew G. Knepley     /* Create initial guess */
3719566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
3729566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
3739566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
3749566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]));
3759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
376adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
377a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
378a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
379a56d3bf6SMatthew G. Knepley 
380a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
3819566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
3829566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
3839566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
3849566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes)));
3859566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r*ce->Nf+f]));
3869566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]));
387adfa7136SMatthew G. Knepley     }
3886af0ca60SMatthew G. Knepley     /* Monitor */
3899566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
39073269098SMatthew G. Knepley     if (!r) {
39173269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
39273269098SMatthew G. Knepley       KSP ksp;
39373269098SMatthew G. Knepley       PC  pc;
39473269098SMatthew G. Knepley 
3959566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
3969566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
3979566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
39873269098SMatthew G. Knepley     }
3996af0ca60SMatthew G. Knepley     /* Cleanup */
4009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
4019566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
4026af0ca60SMatthew G. Knepley   }
4036af0ca60SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
4049566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm[r]));
4056af0ca60SMatthew G. Knepley   }
4066af0ca60SMatthew G. Knepley   /* Fit convergence rate */
4079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr+1, &x, Nr+1, &y));
40846079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
4096af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
4107809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r*ce->Nf+f]);
41146079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
4126af0ca60SMatthew G. Knepley     }
4139566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr+1, x, y, &slope, &intercept));
4146af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
41546079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
41646079b62SMatthew G. Knepley   }
4179566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
4189566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
4192cae373cSMatthew G. Knepley   /* Restore solver */
4209566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
42173269098SMatthew G. Knepley   {
42273269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
42373269098SMatthew G. Knepley     KSP ksp;
42473269098SMatthew G. Knepley     PC  pc;
42573269098SMatthew G. Knepley 
4269566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4279566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4289566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
4299566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
43073269098SMatthew G. Knepley   }
4319566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, ce->idm));
4329566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx));
4339566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4349566063dSJacob Faibussowitsch   PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
435900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
436900f6b5bSMatthew G. Knepley }
437900f6b5bSMatthew G. Knepley 
438900f6b5bSMatthew G. Knepley /*@
439900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
440900f6b5bSMatthew G. Knepley 
441900f6b5bSMatthew G. Knepley   Not collective
442900f6b5bSMatthew G. Knepley 
443900f6b5bSMatthew G. Knepley   Input Parameter:
444900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
445900f6b5bSMatthew G. Knepley 
446900f6b5bSMatthew G. Knepley   Output Parameter:
447900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
448900f6b5bSMatthew G. Knepley 
449900f6b5bSMatthew G. Knepley   Note: The convergence rate alpha is defined by
450900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
451900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
452900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution.
453900f6b5bSMatthew G. Knepley 
454900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
455900f6b5bSMatthew G. Knepley based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
456900f6b5bSMatthew G. Knepley 
457900f6b5bSMatthew G. Knepley   Options database keys:
45867b8a455SSatish Balay + -snes_convergence_estimate - Execute convergence estimation inside SNESSolve() and print out the rate
45967b8a455SSatish Balay - -ts_convergence_estimate - Execute convergence estimation inside TSSolve() and print out the rate
460900f6b5bSMatthew G. Knepley 
461900f6b5bSMatthew G. Knepley   Level: intermediate
462900f6b5bSMatthew G. Knepley 
463900f6b5bSMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
464900f6b5bSMatthew G. Knepley @*/
465900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
466900f6b5bSMatthew G. Knepley {
467900f6b5bSMatthew G. Knepley   PetscInt       f;
468900f6b5bSMatthew G. Knepley 
469900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
4709566063dSJacob Faibussowitsch   if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event));
471900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
4729566063dSJacob Faibussowitsch   PetscCall((*ce->ops->getconvrate)(ce, alpha));
4736af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4746af0ca60SMatthew G. Knepley }
4756af0ca60SMatthew G. Knepley 
4766af0ca60SMatthew G. Knepley /*@
4776af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4786af0ca60SMatthew G. Knepley 
4796af0ca60SMatthew G. Knepley    Collective on SNES
4806af0ca60SMatthew G. Knepley 
4816af0ca60SMatthew G. Knepley    Parameter:
4826af0ca60SMatthew G. Knepley +  snes - iterative context obtained from SNESCreate()
48346079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4846af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4856af0ca60SMatthew G. Knepley 
4866af0ca60SMatthew G. Knepley    Options Database Keys:
4876af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4886af0ca60SMatthew G. Knepley 
4896af0ca60SMatthew G. Knepley    Level: developer
4906af0ca60SMatthew G. Knepley 
4916af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate()
4926af0ca60SMatthew G. Knepley @*/
493a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
4946af0ca60SMatthew G. Knepley {
4956af0ca60SMatthew G. Knepley   PetscBool      isAscii;
4966af0ca60SMatthew G. Knepley 
4976af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii));
4996af0ca60SMatthew G. Knepley   if (isAscii) {
500900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
501a56d3bf6SMatthew G. Knepley 
5029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel));
5039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: "));
5049566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "["));
505a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5069566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
5079566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]));
50846079b62SMatthew G. Knepley     }
5099566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]"));
5109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
5119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel));
5126af0ca60SMatthew G. Knepley   }
5136af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
5146af0ca60SMatthew G. Knepley }
515900f6b5bSMatthew G. Knepley 
516900f6b5bSMatthew G. Knepley /*@
517900f6b5bSMatthew G. Knepley   PetscConvEstCreate - Create a PetscConvEst object
518900f6b5bSMatthew G. Knepley 
519900f6b5bSMatthew G. Knepley   Collective
520900f6b5bSMatthew G. Knepley 
521900f6b5bSMatthew G. Knepley   Input Parameter:
522900f6b5bSMatthew G. Knepley . comm - The communicator for the PetscConvEst object
523900f6b5bSMatthew G. Knepley 
524900f6b5bSMatthew G. Knepley   Output Parameter:
525900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
526900f6b5bSMatthew G. Knepley 
527900f6b5bSMatthew G. Knepley   Level: beginner
528900f6b5bSMatthew G. Knepley 
529900f6b5bSMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
530900f6b5bSMatthew G. Knepley @*/
531900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
532900f6b5bSMatthew G. Knepley {
533900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
534900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
5359566063dSJacob Faibussowitsch   PetscCall(PetscSysInitializePackage());
5369566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView));
537900f6b5bSMatthew G. Knepley   (*ce)->monitor = PETSC_FALSE;
5382e61be88SMatthew G. Knepley   (*ce)->r       = 2.0;
539900f6b5bSMatthew G. Knepley   (*ce)->Nr      = 4;
540900f6b5bSMatthew G. Knepley   (*ce)->event   = -1;
541900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
542900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
543900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
544900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
545900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
546900f6b5bSMatthew G. Knepley }
547