xref: /petsc/src/snes/utils/convest.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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 
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
8d71ae5a4SJacob Faibussowitsch {
96af0ca60SMatthew G. Knepley   PetscInt c;
106af0ca60SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
113ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
126af0ca60SMatthew G. Knepley }
136af0ca60SMatthew G. Knepley 
146af0ca60SMatthew G. Knepley /*@
15420bcc1bSBarry Smith   PetscConvEstDestroy - Destroys a PETSc convergence estimator `PetscConvEst` object
166af0ca60SMatthew G. Knepley 
17c3339decSBarry Smith   Collective
186af0ca60SMatthew G. Knepley 
196af0ca60SMatthew G. Knepley   Input Parameter:
2020f4b53cSBarry Smith . ce - The `PetscConvEst` object
216af0ca60SMatthew G. Knepley 
226af0ca60SMatthew G. Knepley   Level: beginner
236af0ca60SMatthew G. Knepley 
2420f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
256af0ca60SMatthew G. Knepley @*/
26d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
27d71ae5a4SJacob Faibussowitsch {
286af0ca60SMatthew G. Knepley   PetscFunctionBegin;
293ba16761SJacob Faibussowitsch   if (!*ce) PetscFunctionReturn(PETSC_SUCCESS);
30*f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*ce, PETSC_OBJECT_CLASSID, 1);
31*f4f49eeaSPierre Jolivet   if (--((PetscObject)*ce)->refct > 0) {
326af0ca60SMatthew G. Knepley     *ce = NULL;
333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
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));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
396af0ca60SMatthew G. Knepley }
406af0ca60SMatthew G. Knepley 
416af0ca60SMatthew G. Knepley /*@
42420bcc1bSBarry Smith   PetscConvEstSetFromOptions - Sets a convergence estimator `PetscConvEst` object based on values in the options database
436af0ca60SMatthew G. Knepley 
44c3339decSBarry Smith   Collective
456af0ca60SMatthew G. Knepley 
462fe279fdSBarry Smith   Input Parameter:
4720f4b53cSBarry Smith . ce - The `PetscConvEst` object
486af0ca60SMatthew G. Knepley 
496af0ca60SMatthew G. Knepley   Level: beginner
506af0ca60SMatthew G. Knepley 
5120f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
526af0ca60SMatthew G. Knepley @*/
53d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
54d71ae5a4SJacob Faibussowitsch {
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();
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636af0ca60SMatthew G. Knepley }
646af0ca60SMatthew G. Knepley 
656af0ca60SMatthew G. Knepley /*@
6620f4b53cSBarry Smith   PetscConvEstView - Views a `PetscConvEst` object
676af0ca60SMatthew G. Knepley 
68c3339decSBarry Smith   Collective
696af0ca60SMatthew G. Knepley 
706af0ca60SMatthew G. Knepley   Input Parameters:
7120f4b53cSBarry Smith + ce     - The `PetscConvEst` object
72420bcc1bSBarry Smith - viewer - The `PetscViewer`
736af0ca60SMatthew G. Knepley 
746af0ca60SMatthew G. Knepley   Level: beginner
756af0ca60SMatthew G. Knepley 
7620f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscViewer`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
776af0ca60SMatthew G. Knepley @*/
78d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
79d71ae5a4SJacob Faibussowitsch {
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));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
846af0ca60SMatthew G. Knepley }
856af0ca60SMatthew G. Knepley 
866af0ca60SMatthew G. Knepley /*@
876af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
886af0ca60SMatthew G. Knepley 
8920f4b53cSBarry Smith   Not Collective
906af0ca60SMatthew G. Knepley 
916af0ca60SMatthew G. Knepley   Input Parameter:
9220f4b53cSBarry Smith . 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 
9920f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1006af0ca60SMatthew G. Knepley @*/
101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
102d71ae5a4SJacob Faibussowitsch {
1036af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1046af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
1054f572ea9SToby Isaac   PetscAssertPointer(solver, 2);
106900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1086af0ca60SMatthew G. Knepley }
1096af0ca60SMatthew G. Knepley 
1106af0ca60SMatthew G. Knepley /*@
1116af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1126af0ca60SMatthew G. Knepley 
11320f4b53cSBarry Smith   Not Collective
1146af0ca60SMatthew G. Knepley 
1156af0ca60SMatthew G. Knepley   Input Parameters:
11620f4b53cSBarry Smith + ce     - The `PetscConvEst` object
117420bcc1bSBarry Smith - solver - The solver, must be a `KSP`, `SNES`, or `TS` object with an attached `DM`/`DS`, that can compute an exact solution
1186af0ca60SMatthew G. Knepley 
1196af0ca60SMatthew G. Knepley   Level: intermediate
1206af0ca60SMatthew G. Knepley 
12120f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1226af0ca60SMatthew G. Knepley @*/
123d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
124d71ae5a4SJacob Faibussowitsch {
1256af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1266af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
127900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
128900f6b5bSMatthew G. Knepley   ce->solver = solver;
129dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, setsolver, solver);
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1316af0ca60SMatthew G. Knepley }
1326af0ca60SMatthew G. Knepley 
1336af0ca60SMatthew G. Knepley /*@
134420bcc1bSBarry Smith   PetscConvEstSetUp - After the solver is specified, create data structures needed for estimating convergence
1356af0ca60SMatthew G. Knepley 
136c3339decSBarry Smith   Collective
1376af0ca60SMatthew G. Knepley 
1382fe279fdSBarry Smith   Input Parameter:
13920f4b53cSBarry Smith . ce - The `PetscConvEst` object
1406af0ca60SMatthew G. Knepley 
1416af0ca60SMatthew G. Knepley   Level: beginner
1426af0ca60SMatthew G. Knepley 
14320f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1446af0ca60SMatthew G. Knepley @*/
145d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
146d71ae5a4SJacob Faibussowitsch {
147083401c6SMatthew G. Knepley   PetscInt Nf, f, Nds, s;
1486af0ca60SMatthew G. Knepley 
1496af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1509566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(ce->idm, &Nf));
151900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2((ce->Nr + 1) * ce->Nf, &ce->dofs, (ce->Nr + 1) * ce->Nf, &ce->errors));
1539566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs));
154900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
1559566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(ce->idm, &Nds));
156083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
157083401c6SMatthew G. Knepley     PetscDS         ds;
158083401c6SMatthew G. Knepley     DMLabel         label;
159083401c6SMatthew G. Knepley     IS              fieldIS;
160083401c6SMatthew G. Knepley     const PetscInt *fields;
161083401c6SMatthew G. Knepley     PetscInt        dsNf;
162083401c6SMatthew G. Knepley 
16307218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds, NULL));
1649566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &dsNf));
1659566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
166083401c6SMatthew G. Knepley     for (f = 0; f < dsNf; ++f) {
167083401c6SMatthew G. Knepley       const PetscInt field = fields[f];
1689566063dSJacob Faibussowitsch       PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]));
169083401c6SMatthew G. Knepley     }
1709566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
171083401c6SMatthew G. Knepley   }
172ad540459SPierre Jolivet   for (f = 0; f < Nf; ++f) 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);
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1746af0ca60SMatthew G. Knepley }
1756af0ca60SMatthew G. Knepley 
176d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
177d71ae5a4SJacob Faibussowitsch {
1786af0ca60SMatthew G. Knepley   PetscFunctionBegin;
179900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
180900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
181900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
182dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, initguess, r, dm, u);
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184900f6b5bSMatthew G. Knepley }
185900f6b5bSMatthew G. Knepley 
186d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
187d71ae5a4SJacob Faibussowitsch {
188900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
189900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
190900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
191900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1924f572ea9SToby Isaac   PetscAssertPointer(errors, 5);
193dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, computeerror, r, dm, u, errors);
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
195900f6b5bSMatthew G. Knepley }
196900f6b5bSMatthew G. Knepley 
197f2cacb80SMatthew G. Knepley /*@
198f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
199f2cacb80SMatthew G. Knepley 
200c3339decSBarry Smith   Collective
201f2cacb80SMatthew G. Knepley 
202d8d19677SJose E. Roman   Input Parameters:
203f6dfbefdSBarry Smith + ce - The `PetscConvEst` object
204f2cacb80SMatthew G. Knepley - r  - The refinement level
205f2cacb80SMatthew G. Knepley 
206420bcc1bSBarry Smith   Options Database Key:
207ee300463SSatish Balay . -convest_monitor - Activate the monitor
208f2cacb80SMatthew G. Knepley 
209f2cacb80SMatthew G. Knepley   Level: intermediate
210f2cacb80SMatthew G. Knepley 
211f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
212f2cacb80SMatthew G. Knepley @*/
213d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
214d71ae5a4SJacob Faibussowitsch {
215900f6b5bSMatthew G. Knepley   MPI_Comm comm;
216900f6b5bSMatthew G. Knepley   PetscInt f;
217900f6b5bSMatthew G. Knepley 
218900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
219900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2207809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r * ce->Nf];
221900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r * ce->Nf];
222900f6b5bSMatthew G. Knepley 
2239566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)ce, &comm));
2249566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "N: "));
2259566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
2267809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2279566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
22863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f]));
2297809adefSMatthew G. Knepley     }
2309566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2319566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "  "));
2329566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Error: "));
2339566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
234900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2359566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2369566063dSJacob Faibussowitsch       if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11"));
2379566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "%g", (double)errors[f]));
238900f6b5bSMatthew G. Knepley     }
2399566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n"));
241900f6b5bSMatthew G. Knepley   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243900f6b5bSMatthew G. Knepley }
244900f6b5bSMatthew G. Knepley 
245d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
246d71ae5a4SJacob Faibussowitsch {
247900f6b5bSMatthew G. Knepley   PetscClassId id;
248900f6b5bSMatthew G. Knepley 
249900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
25108401ef6SPierre Jolivet   PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
2529566063dSJacob Faibussowitsch   PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm));
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254900f6b5bSMatthew G. Knepley }
255900f6b5bSMatthew G. Knepley 
256d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
257d71ae5a4SJacob Faibussowitsch {
258900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u));
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
261900f6b5bSMatthew G. Knepley }
262900f6b5bSMatthew G. Knepley 
263d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
264d71ae5a4SJacob Faibussowitsch {
265900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors));
2673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
268900f6b5bSMatthew G. Knepley }
269900f6b5bSMatthew G. Knepley 
2708c401167SJed Brown static PetscErrorCode PetscConvEstSetJacobianNullSpace_Private(PetscConvEst ce, SNES snes)
271d71ae5a4SJacob Faibussowitsch {
272478db826SMatthew G. Knepley   DM       dm;
273478db826SMatthew G. Knepley   PetscInt f;
274478db826SMatthew G. Knepley 
275478db826SMatthew G. Knepley   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
277478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
278478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
279478db826SMatthew G. Knepley 
2809566063dSJacob Faibussowitsch     PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr));
281478db826SMatthew G. Knepley     if (nspconstr) {
282478db826SMatthew G. Knepley       MatNullSpace nullsp;
283478db826SMatthew G. Knepley       Mat          J;
284478db826SMatthew G. Knepley 
2859566063dSJacob Faibussowitsch       PetscCall((*nspconstr)(dm, f, f, &nullsp));
2869566063dSJacob Faibussowitsch       PetscCall(SNESSetUp(snes));
2879566063dSJacob Faibussowitsch       PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
2889566063dSJacob Faibussowitsch       PetscCall(MatSetNullSpace(J, nullsp));
2899566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullsp));
290478db826SMatthew G. Knepley       break;
291478db826SMatthew G. Knepley     }
292478db826SMatthew G. Knepley   }
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294478db826SMatthew G. Knepley }
295478db826SMatthew G. Knepley 
296d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
297d71ae5a4SJacob Faibussowitsch {
298900f6b5bSMatthew G. Knepley   SNES        snes = (SNES)ce->solver;
299900f6b5bSMatthew G. Knepley   DM         *dm;
300900f6b5bSMatthew G. Knepley   PetscObject disc;
301900f6b5bSMatthew G. Knepley   PetscReal  *x, *y, slope, intercept;
3027809adefSMatthew G. Knepley   PetscInt    Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
303900f6b5bSMatthew G. Knepley   void       *ctx;
304900f6b5bSMatthew G. Knepley 
305900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
30608401ef6SPierre Jolivet   PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r);
3079566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
3089566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
3099566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
3109566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
3119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr + 1), &dm));
3126af0ca60SMatthew G. Knepley   /* Loop over meshes */
313900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3146af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
315900f6b5bSMatthew G. Knepley     Vec           u;
316adfa7136SMatthew G. Knepley     PetscLogStage stage;
317adfa7136SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
318900f6b5bSMatthew G. Knepley     const char   *dmname, *uname;
319adfa7136SMatthew G. Knepley 
32063a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r));
3219566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
3229566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
3239566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3246af0ca60SMatthew G. Knepley     if (r > 0) {
325b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
3269566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r]));
3279566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r - 1]));
328b2df8587SMatthew G. Knepley       } else {
329b2df8587SMatthew G. Knepley         DM cdm, rcdm;
330b2df8587SMatthew G. Knepley 
3319566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r - 1], &dm[r]));
3329566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r - 1], dm[r]));
3339566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm));
3349566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r], &rcdm));
3359566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
336b2df8587SMatthew G. Knepley       }
3379566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
3389566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname));
3399566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname));
340478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3418cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
342900f6b5bSMatthew G. Knepley 
3439566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr));
3449566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr));
3458f892730SMatthew G. Knepley       }
3466af0ca60SMatthew G. Knepley     }
3479566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
3486af0ca60SMatthew G. Knepley     /* Create solution */
3499566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
3509566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
3519566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
3529566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)u, uname));
3536af0ca60SMatthew G. Knepley     /* Setup solver */
3549566063dSJacob Faibussowitsch     PetscCall(SNESReset(snes));
3559566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm[r]));
3566493148fSStefano Zampini     PetscCall(DMPlexSetSNESLocalFEM(dm[r], PETSC_FALSE, ctx));
3579566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
358478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
3598c401167SJed Brown     PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes));
3606af0ca60SMatthew G. Knepley     /* Create initial guess */
3619566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
3629566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
3639566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
3649566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf]));
3659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
366adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
367a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
368a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
369a56d3bf6SMatthew G. Knepley 
370a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
3719566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
3729566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
3739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
374712fec58SPierre Jolivet       PetscCall(MPIU_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
3759566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f]));
3769566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f]));
377adfa7136SMatthew G. Knepley     }
3786af0ca60SMatthew G. Knepley     /* Monitor */
3799566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
38073269098SMatthew G. Knepley     if (!r) {
38173269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
38273269098SMatthew G. Knepley       KSP ksp;
38373269098SMatthew G. Knepley       PC  pc;
38473269098SMatthew G. Knepley 
3859566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
3869566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
3879566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
38873269098SMatthew G. Knepley     }
3896af0ca60SMatthew G. Knepley     /* Cleanup */
3909566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
3919566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
3926af0ca60SMatthew G. Knepley   }
39348a46eb9SPierre Jolivet   for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r]));
3946af0ca60SMatthew G. Knepley   /* Fit convergence rate */
3959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y));
39646079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
3976af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
3987809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]);
39946079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]);
4006af0ca60SMatthew G. Knepley     }
4019566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept));
4026af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
40346079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
40446079b62SMatthew G. Knepley   }
4059566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
4072cae373cSMatthew G. Knepley   /* Restore solver */
4089566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
40973269098SMatthew G. Knepley   {
41073269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
41173269098SMatthew G. Knepley     KSP ksp;
41273269098SMatthew G. Knepley     PC  pc;
41373269098SMatthew G. Knepley 
4149566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4159566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4169566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
4179566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
41873269098SMatthew G. Knepley   }
4199566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, ce->idm));
4206493148fSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(ce->idm, PETSC_FALSE, ctx));
4219566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4228c401167SJed Brown   PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes));
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424900f6b5bSMatthew G. Knepley }
425900f6b5bSMatthew G. Knepley 
426900f6b5bSMatthew G. Knepley /*@
427900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
428900f6b5bSMatthew G. Knepley 
42920f4b53cSBarry Smith   Not Collective
430900f6b5bSMatthew G. Knepley 
431900f6b5bSMatthew G. Knepley   Input Parameter:
432f6dfbefdSBarry Smith . ce - The `PetscConvEst` object
433900f6b5bSMatthew G. Knepley 
434900f6b5bSMatthew G. Knepley   Output Parameter:
435900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
436900f6b5bSMatthew G. Knepley 
43720f4b53cSBarry Smith   Options Database Keys:
438f6dfbefdSBarry Smith + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate
439f6dfbefdSBarry Smith - -ts_convergence_estimate   - Execute convergence estimation inside `TSSolve()` and print out the rate
440f6dfbefdSBarry Smith 
44120f4b53cSBarry Smith   Level: intermediate
44220f4b53cSBarry Smith 
443f6dfbefdSBarry Smith   Notes:
444f6dfbefdSBarry Smith   The convergence rate alpha is defined by
445900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
446900f6b5bSMatthew G. Knepley   where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
447900f6b5bSMatthew G. Knepley   spatial resolution and \Delta t for the temporal resolution.
448900f6b5bSMatthew G. Knepley 
449900f6b5bSMatthew G. Knepley   We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
45020f4b53cSBarry Smith   based upon the exact solution in the `PetscDS`, and then fit the result to our model above using linear regression.
451900f6b5bSMatthew G. Knepley 
45242747ad1SJacob Faibussowitsch .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `SNESSolve()`, `TSSolve()`
453900f6b5bSMatthew G. Knepley @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
455d71ae5a4SJacob Faibussowitsch {
456900f6b5bSMatthew G. Knepley   PetscInt f;
457900f6b5bSMatthew G. Knepley 
458900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event));
460900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
461dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, getconvrate, alpha);
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4636af0ca60SMatthew G. Knepley }
4646af0ca60SMatthew G. Knepley 
4656af0ca60SMatthew G. Knepley /*@
466420bcc1bSBarry Smith   PetscConvEstRateView - Displays the convergence rate obtained from `PetscConvEstGetConvRate()` using a `PetscViewer`
4676af0ca60SMatthew G. Knepley 
468c3339decSBarry Smith   Collective
4696af0ca60SMatthew G. Knepley 
470e4094ef1SJacob Faibussowitsch   Input Parameters:
471e4094ef1SJacob Faibussowitsch + ce     - iterative context obtained from `SNESCreate()`
47246079b62SMatthew G. Knepley . alpha  - the convergence rate for each field
4736af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason
4746af0ca60SMatthew G. Knepley 
47520f4b53cSBarry Smith   Options Database Key:
4766af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate
4776af0ca60SMatthew G. Knepley 
4786af0ca60SMatthew G. Knepley   Level: developer
4796af0ca60SMatthew G. Knepley 
480420bcc1bSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetConvRate()`
4816af0ca60SMatthew G. Knepley @*/
482d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
483d71ae5a4SJacob Faibussowitsch {
4846af0ca60SMatthew G. Knepley   PetscBool isAscii;
4856af0ca60SMatthew G. Knepley 
4866af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
4886af0ca60SMatthew G. Knepley   if (isAscii) {
489900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
490a56d3bf6SMatthew G. Knepley 
4919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel));
4929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: "));
4939566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "["));
494a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
4959566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
4969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f]));
49746079b62SMatthew G. Knepley     }
4989566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]"));
4999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
5009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel));
5016af0ca60SMatthew G. Knepley   }
5023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5036af0ca60SMatthew G. Knepley }
504900f6b5bSMatthew G. Knepley 
505900f6b5bSMatthew G. Knepley /*@
506f6dfbefdSBarry Smith   PetscConvEstCreate - Create a `PetscConvEst` object
507900f6b5bSMatthew G. Knepley 
508900f6b5bSMatthew G. Knepley   Collective
509900f6b5bSMatthew G. Knepley 
510900f6b5bSMatthew G. Knepley   Input Parameter:
511f6dfbefdSBarry Smith . comm - The communicator for the `PetscConvEst` object
512900f6b5bSMatthew G. Knepley 
513900f6b5bSMatthew G. Knepley   Output Parameter:
514f6dfbefdSBarry Smith . ce - The `PetscConvEst` object
515900f6b5bSMatthew G. Knepley 
516900f6b5bSMatthew G. Knepley   Level: beginner
517900f6b5bSMatthew G. Knepley 
518420bcc1bSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()`, `DMAdaptorCreate()`, `DMAdaptor`
519900f6b5bSMatthew G. Knepley @*/
520d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
521d71ae5a4SJacob Faibussowitsch {
522900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
5234f572ea9SToby Isaac   PetscAssertPointer(ce, 2);
5249566063dSJacob Faibussowitsch   PetscCall(PetscSysInitializePackage());
5259566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView));
526900f6b5bSMatthew G. Knepley   (*ce)->monitor           = PETSC_FALSE;
5272e61be88SMatthew G. Knepley   (*ce)->r                 = 2.0;
528900f6b5bSMatthew G. Knepley   (*ce)->Nr                = 4;
529900f6b5bSMatthew G. Knepley   (*ce)->event             = -1;
530900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
531900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
532900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
533900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535900f6b5bSMatthew G. Knepley }
536