xref: /petsc/src/snes/utils/convest.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
24db781477SPatrick Sanan .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 
51db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
526af0ca60SMatthew G. Knepley @*/
536af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
546af0ca60SMatthew G. Knepley {
556af0ca60SMatthew G. Knepley   PetscFunctionBegin;
56d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL));
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL));
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL));
61d0609cedSBarry Smith   PetscOptionsEnd();
626af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
636af0ca60SMatthew G. Knepley }
646af0ca60SMatthew G. Knepley 
656af0ca60SMatthew G. Knepley /*@
666af0ca60SMatthew G. Knepley   PetscConvEstView - Views a PetscConvEst object
676af0ca60SMatthew G. Knepley 
686af0ca60SMatthew G. Knepley   Collective on PetscConvEst
696af0ca60SMatthew G. Knepley 
706af0ca60SMatthew G. Knepley   Input Parameters:
716af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
726af0ca60SMatthew G. Knepley - viewer - The PetscViewer object
736af0ca60SMatthew G. Knepley 
746af0ca60SMatthew G. Knepley   Level: beginner
756af0ca60SMatthew G. Knepley 
76db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
776af0ca60SMatthew G. Knepley @*/
786af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
796af0ca60SMatthew G. Knepley {
806af0ca60SMatthew G. Knepley   PetscFunctionBegin;
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer));
8263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %" PetscInt_FMT " levels\n", ce->Nr+1));
836af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
846af0ca60SMatthew G. Knepley }
856af0ca60SMatthew G. Knepley 
866af0ca60SMatthew G. Knepley /*@
876af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
886af0ca60SMatthew G. Knepley 
896af0ca60SMatthew G. Knepley   Not collective
906af0ca60SMatthew G. Knepley 
916af0ca60SMatthew G. Knepley   Input Parameter:
926af0ca60SMatthew G. Knepley . ce     - The PetscConvEst object
936af0ca60SMatthew G. Knepley 
946af0ca60SMatthew G. Knepley   Output Parameter:
95900f6b5bSMatthew G. Knepley . solver - The solver
966af0ca60SMatthew G. Knepley 
976af0ca60SMatthew G. Knepley   Level: intermediate
986af0ca60SMatthew G. Knepley 
99db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1006af0ca60SMatthew G. Knepley @*/
101900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
1026af0ca60SMatthew G. Knepley {
1036af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1046af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
105900f6b5bSMatthew G. Knepley   PetscValidPointer(solver, 2);
106900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1076af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1086af0ca60SMatthew G. Knepley }
1096af0ca60SMatthew G. Knepley 
1106af0ca60SMatthew G. Knepley /*@
1116af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1126af0ca60SMatthew G. Knepley 
1136af0ca60SMatthew G. Knepley   Not collective
1146af0ca60SMatthew G. Knepley 
1156af0ca60SMatthew G. Knepley   Input Parameters:
1166af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
117900f6b5bSMatthew G. Knepley - solver - The solver
1186af0ca60SMatthew G. Knepley 
1196af0ca60SMatthew G. Knepley   Level: intermediate
1206af0ca60SMatthew G. Knepley 
1216af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1226af0ca60SMatthew G. Knepley 
123db781477SPatrick Sanan .seealso: `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1246af0ca60SMatthew G. Knepley @*/
125900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
1266af0ca60SMatthew G. Knepley {
1276af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1286af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
129900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
130900f6b5bSMatthew G. Knepley   ce->solver = solver;
131*dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce,setsolver , solver);
1326af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1336af0ca60SMatthew G. Knepley }
1346af0ca60SMatthew G. Knepley 
1356af0ca60SMatthew G. Knepley /*@
1366af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1376af0ca60SMatthew G. Knepley 
1386af0ca60SMatthew G. Knepley   Collective on PetscConvEst
1396af0ca60SMatthew G. Knepley 
1406af0ca60SMatthew G. Knepley   Input Parameters:
1416af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1426af0ca60SMatthew G. Knepley 
1436af0ca60SMatthew G. Knepley   Level: beginner
1446af0ca60SMatthew G. Knepley 
145db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1466af0ca60SMatthew G. Knepley @*/
1470955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
1486af0ca60SMatthew G. Knepley {
149083401c6SMatthew G. Knepley   PetscInt       Nf, f, Nds, s;
1506af0ca60SMatthew G. Knepley 
1516af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(ce->idm, &Nf));
153900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2((ce->Nr+1)*ce->Nf, &ce->dofs, (ce->Nr+1)*ce->Nf, &ce->errors));
1559566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs));
156900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
1579566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(ce->idm, &Nds));
158083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
159083401c6SMatthew G. Knepley     PetscDS         ds;
160083401c6SMatthew G. Knepley     DMLabel         label;
161083401c6SMatthew G. Knepley     IS              fieldIS;
162083401c6SMatthew G. Knepley     const PetscInt *fields;
163083401c6SMatthew G. Knepley     PetscInt        dsNf;
164083401c6SMatthew G. Knepley 
1659566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds));
1669566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &dsNf));
1679566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
168083401c6SMatthew G. Knepley     for (f = 0; f < dsNf; ++f) {
169083401c6SMatthew G. Knepley       const PetscInt field = fields[f];
1709566063dSJacob Faibussowitsch       PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]));
171083401c6SMatthew G. Knepley     }
1729566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
173083401c6SMatthew G. Knepley   }
174900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
17563a3b9bcSJacob Faibussowitsch     PetscCheck(ce->exactSol[f],PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %" PetscInt_FMT, f);
1766af0ca60SMatthew G. Knepley   }
1776af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1786af0ca60SMatthew G. Knepley }
1796af0ca60SMatthew G. Knepley 
180900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
1816af0ca60SMatthew G. Knepley {
1826af0ca60SMatthew G. Knepley   PetscFunctionBegin;
183900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
184900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
185900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
186*dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce,initguess , r, dm, u);
187900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
188900f6b5bSMatthew G. Knepley }
189900f6b5bSMatthew G. Knepley 
190900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
191900f6b5bSMatthew G. Knepley {
192900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
193900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
194900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
195900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
196900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
197*dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce,computeerror , r, dm, u, errors);
198900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
199900f6b5bSMatthew G. Knepley }
200900f6b5bSMatthew G. Knepley 
201f2cacb80SMatthew G. Knepley /*@
202f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
203f2cacb80SMatthew G. Knepley 
204f2cacb80SMatthew G. Knepley   Collective on PetscConvEst
205f2cacb80SMatthew G. Knepley 
206d8d19677SJose E. Roman   Input Parameters:
207f2cacb80SMatthew G. Knepley + ce - The PetscConvEst object
208f2cacb80SMatthew G. Knepley - r  - The refinement level
209f2cacb80SMatthew G. Knepley 
210f2cacb80SMatthew G. Knepley   Options database keys:
211ee300463SSatish Balay . -convest_monitor - Activate the monitor
212f2cacb80SMatthew G. Knepley 
213f2cacb80SMatthew G. Knepley   Level: intermediate
214f2cacb80SMatthew G. Knepley 
215db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
216f2cacb80SMatthew G. Knepley @*/
217f2cacb80SMatthew G. Knepley PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
218900f6b5bSMatthew G. Knepley {
219900f6b5bSMatthew G. Knepley   MPI_Comm       comm;
220900f6b5bSMatthew G. Knepley   PetscInt       f;
221900f6b5bSMatthew G. Knepley 
222900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
223900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2247809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r*ce->Nf];
225900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r*ce->Nf];
226900f6b5bSMatthew G. Knepley 
2279566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject) ce, &comm));
2289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "N: "));
2299566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
2307809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2319566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
23263a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f]));
2337809adefSMatthew G. Knepley     }
2349566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2359566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "  "));
2369566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Error: "));
2379566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
238900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2399566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2409566063dSJacob Faibussowitsch       if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11"));
2419566063dSJacob Faibussowitsch       else                     PetscCall(PetscPrintf(comm, "%g", (double) errors[f]));
242900f6b5bSMatthew G. Knepley     }
2439566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n"));
245900f6b5bSMatthew G. Knepley   }
246900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
247900f6b5bSMatthew G. Knepley }
248900f6b5bSMatthew G. Knepley 
249900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
250900f6b5bSMatthew G. Knepley {
251900f6b5bSMatthew G. Knepley   PetscClassId   id;
252900f6b5bSMatthew G. Knepley 
253900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2549566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
25508401ef6SPierre Jolivet   PetscCheck(id == SNES_CLASSID,PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
2569566063dSJacob Faibussowitsch   PetscCall(SNESGetDM((SNES) ce->solver, &ce->idm));
257900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
258900f6b5bSMatthew G. Knepley }
259900f6b5bSMatthew G. Knepley 
260900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
261900f6b5bSMatthew G. Knepley {
262900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u));
264900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
265900f6b5bSMatthew G. Knepley }
266900f6b5bSMatthew G. Knepley 
267900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
268900f6b5bSMatthew G. Knepley {
269900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2709566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors));
271900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
272900f6b5bSMatthew G. Knepley }
273900f6b5bSMatthew G. Knepley 
274478db826SMatthew G. Knepley static PetscErrorCode PetscConvEstSetJacobianNullspace_Private(PetscConvEst ce, SNES snes)
275478db826SMatthew G. Knepley {
276478db826SMatthew G. Knepley   DM             dm;
277478db826SMatthew G. Knepley   PetscInt       f;
278478db826SMatthew G. Knepley 
279478db826SMatthew G. Knepley   PetscFunctionBegin;
2809566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
281478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
282478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
283478db826SMatthew G. Knepley 
2849566063dSJacob Faibussowitsch     PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr));
285478db826SMatthew G. Knepley     if (nspconstr) {
286478db826SMatthew G. Knepley       MatNullSpace nullsp;
287478db826SMatthew G. Knepley       Mat          J;
288478db826SMatthew G. Knepley 
2899566063dSJacob Faibussowitsch       PetscCall((*nspconstr)(dm, f, f,&nullsp));
2909566063dSJacob Faibussowitsch       PetscCall(SNESSetUp(snes));
2919566063dSJacob Faibussowitsch       PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
2929566063dSJacob Faibussowitsch       PetscCall(MatSetNullSpace(J, nullsp));
2939566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullsp));
294478db826SMatthew G. Knepley       break;
295478db826SMatthew G. Knepley     }
296478db826SMatthew G. Knepley   }
297478db826SMatthew G. Knepley   PetscFunctionReturn(0);
298478db826SMatthew G. Knepley }
299478db826SMatthew G. Knepley 
300900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
301900f6b5bSMatthew G. Knepley {
302900f6b5bSMatthew G. Knepley   SNES           snes = (SNES) ce->solver;
303900f6b5bSMatthew G. Knepley   DM            *dm;
304900f6b5bSMatthew G. Knepley   PetscObject    disc;
305900f6b5bSMatthew G. Knepley   PetscReal     *x, *y, slope, intercept;
3067809adefSMatthew G. Knepley   PetscInt       Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
307900f6b5bSMatthew G. Knepley   void          *ctx;
308900f6b5bSMatthew G. Knepley 
309900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
31008401ef6SPierre Jolivet   PetscCheck(ce->r == 2.0,PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
3119566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
3129566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
3139566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
3149566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
3159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr+1), &dm));
3166af0ca60SMatthew G. Knepley   /* Loop over meshes */
317900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3186af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
319900f6b5bSMatthew G. Knepley     Vec           u;
320e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG)
321adfa7136SMatthew G. Knepley     PetscLogStage stage;
322e5ed2c37SJose E. Roman #endif
323adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
324900f6b5bSMatthew G. Knepley     const char   *dmname, *uname;
325adfa7136SMatthew G. Knepley 
32663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %" PetscInt_FMT, r));
327608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
3289566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
3299566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
330608e5a7aSMatthew G. Knepley #endif
3319566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3326af0ca60SMatthew G. Knepley     if (r > 0) {
333b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
3349566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]));
3359566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r-1]));
336b2df8587SMatthew G. Knepley       } else {
337b2df8587SMatthew G. Knepley         DM cdm, rcdm;
338b2df8587SMatthew G. Knepley 
3399566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r-1], &dm[r]));
3409566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r-1], dm[r]));
3419566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r-1], &cdm));
3429566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r],   &rcdm));
3439566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
344b2df8587SMatthew G. Knepley       }
3459566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
3469566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) dm[r-1], &dmname));
3479566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) dm[r], dmname));
348478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3498cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
350900f6b5bSMatthew G. Knepley 
3519566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr));
3529566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r],   f,  nspconstr));
3538f892730SMatthew G. Knepley       }
3546af0ca60SMatthew G. Knepley     }
3559566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
3566af0ca60SMatthew G. Knepley     /* Create solution */
3579566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
3589566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
3599566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
3609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) u, uname));
3616af0ca60SMatthew G. Knepley     /* Setup solver */
3629566063dSJacob Faibussowitsch     PetscCall(SNESReset(snes));
3639566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm[r]));
3649566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx));
3659566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
366478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
3679566063dSJacob Faibussowitsch     PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
3686af0ca60SMatthew G. Knepley     /* Create initial guess */
3699566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
3709566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
3719566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
3729566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]));
3739566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
374adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
375a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
376a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
377a56d3bf6SMatthew G. Knepley 
378a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
3799566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
3809566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
3819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
3829566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes)));
3839566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r*ce->Nf+f]));
3849566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]));
385adfa7136SMatthew G. Knepley     }
3866af0ca60SMatthew G. Knepley     /* Monitor */
3879566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
38873269098SMatthew G. Knepley     if (!r) {
38973269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
39073269098SMatthew G. Knepley       KSP ksp;
39173269098SMatthew G. Knepley       PC  pc;
39273269098SMatthew G. Knepley 
3939566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
3949566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
3959566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
39673269098SMatthew G. Knepley     }
3976af0ca60SMatthew G. Knepley     /* Cleanup */
3989566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
3999566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
4006af0ca60SMatthew G. Knepley   }
4016af0ca60SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
4029566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm[r]));
4036af0ca60SMatthew G. Knepley   }
4046af0ca60SMatthew G. Knepley   /* Fit convergence rate */
4059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr+1, &x, Nr+1, &y));
40646079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
4076af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
4087809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r*ce->Nf+f]);
40946079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
4106af0ca60SMatthew G. Knepley     }
4119566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr+1, x, y, &slope, &intercept));
4126af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
41346079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
41446079b62SMatthew G. Knepley   }
4159566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
4169566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
4172cae373cSMatthew G. Knepley   /* Restore solver */
4189566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
41973269098SMatthew G. Knepley   {
42073269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
42173269098SMatthew G. Knepley     KSP ksp;
42273269098SMatthew G. Knepley     PC  pc;
42373269098SMatthew G. Knepley 
4249566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4259566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4269566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
4279566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
42873269098SMatthew G. Knepley   }
4299566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, ce->idm));
4309566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx));
4319566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4329566063dSJacob Faibussowitsch   PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
433900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
434900f6b5bSMatthew G. Knepley }
435900f6b5bSMatthew G. Knepley 
436900f6b5bSMatthew G. Knepley /*@
437900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
438900f6b5bSMatthew G. Knepley 
439900f6b5bSMatthew G. Knepley   Not collective
440900f6b5bSMatthew G. Knepley 
441900f6b5bSMatthew G. Knepley   Input Parameter:
442900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
443900f6b5bSMatthew G. Knepley 
444900f6b5bSMatthew G. Knepley   Output Parameter:
445900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
446900f6b5bSMatthew G. Knepley 
447900f6b5bSMatthew G. Knepley   Note: The convergence rate alpha is defined by
448900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
449900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
450900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution.
451900f6b5bSMatthew G. Knepley 
452900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
453900f6b5bSMatthew G. Knepley based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
454900f6b5bSMatthew G. Knepley 
455900f6b5bSMatthew G. Knepley   Options database keys:
45667b8a455SSatish Balay + -snes_convergence_estimate - Execute convergence estimation inside SNESSolve() and print out the rate
45767b8a455SSatish Balay - -ts_convergence_estimate - Execute convergence estimation inside TSSolve() and print out the rate
458900f6b5bSMatthew G. Knepley 
459900f6b5bSMatthew G. Knepley   Level: intermediate
460900f6b5bSMatthew G. Knepley 
461db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
462900f6b5bSMatthew G. Knepley @*/
463900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
464900f6b5bSMatthew G. Knepley {
465900f6b5bSMatthew G. Knepley   PetscInt       f;
466900f6b5bSMatthew G. Knepley 
467900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event));
469900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
470*dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce,getconvrate , alpha);
4716af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4726af0ca60SMatthew G. Knepley }
4736af0ca60SMatthew G. Knepley 
4746af0ca60SMatthew G. Knepley /*@
4756af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4766af0ca60SMatthew G. Knepley 
4776af0ca60SMatthew G. Knepley    Collective on SNES
4786af0ca60SMatthew G. Knepley 
4796af0ca60SMatthew G. Knepley    Parameter:
4806af0ca60SMatthew G. Knepley +  snes - iterative context obtained from SNESCreate()
48146079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4826af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4836af0ca60SMatthew G. Knepley 
4846af0ca60SMatthew G. Knepley    Options Database Keys:
4856af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4866af0ca60SMatthew G. Knepley 
4876af0ca60SMatthew G. Knepley    Level: developer
4886af0ca60SMatthew G. Knepley 
489db781477SPatrick Sanan .seealso: `PetscConvEstGetRate()`
4906af0ca60SMatthew G. Knepley @*/
491a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
4926af0ca60SMatthew G. Knepley {
4936af0ca60SMatthew G. Knepley   PetscBool      isAscii;
4946af0ca60SMatthew G. Knepley 
4956af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii));
4976af0ca60SMatthew G. Knepley   if (isAscii) {
498900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
499a56d3bf6SMatthew G. Knepley 
5009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel));
5019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: "));
5029566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "["));
503a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5049566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
5059566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]));
50646079b62SMatthew G. Knepley     }
5079566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]"));
5089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
5099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel));
5106af0ca60SMatthew G. Knepley   }
5116af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
5126af0ca60SMatthew G. Knepley }
513900f6b5bSMatthew G. Knepley 
514900f6b5bSMatthew G. Knepley /*@
515900f6b5bSMatthew G. Knepley   PetscConvEstCreate - Create a PetscConvEst object
516900f6b5bSMatthew G. Knepley 
517900f6b5bSMatthew G. Knepley   Collective
518900f6b5bSMatthew G. Knepley 
519900f6b5bSMatthew G. Knepley   Input Parameter:
520900f6b5bSMatthew G. Knepley . comm - The communicator for the PetscConvEst object
521900f6b5bSMatthew G. Knepley 
522900f6b5bSMatthew G. Knepley   Output Parameter:
523900f6b5bSMatthew G. Knepley . ce   - The PetscConvEst object
524900f6b5bSMatthew G. Knepley 
525900f6b5bSMatthew G. Knepley   Level: beginner
526900f6b5bSMatthew G. Knepley 
527db781477SPatrick Sanan .seealso: `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()`
528900f6b5bSMatthew G. Knepley @*/
529900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
530900f6b5bSMatthew G. Knepley {
531900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
532900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
5339566063dSJacob Faibussowitsch   PetscCall(PetscSysInitializePackage());
5349566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView));
535900f6b5bSMatthew G. Knepley   (*ce)->monitor = PETSC_FALSE;
5362e61be88SMatthew G. Knepley   (*ce)->r       = 2.0;
537900f6b5bSMatthew G. Knepley   (*ce)->Nr      = 4;
538900f6b5bSMatthew G. Knepley   (*ce)->event   = -1;
539900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
540900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
541900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
542900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
543900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
544900f6b5bSMatthew G. Knepley }
545