xref: /petsc/src/snes/utils/convest.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
79371c9d4SSatish Balay static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
86af0ca60SMatthew G. Knepley   PetscInt c;
96af0ca60SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
106af0ca60SMatthew G. Knepley   return 0;
116af0ca60SMatthew G. Knepley }
126af0ca60SMatthew G. Knepley 
136af0ca60SMatthew G. Knepley /*@
146af0ca60SMatthew G. Knepley   PetscConvEstDestroy - Destroys a PetscConvEst object
156af0ca60SMatthew G. Knepley 
16*f6dfbefdSBarry Smith   Collective on ce
176af0ca60SMatthew G. Knepley 
186af0ca60SMatthew G. Knepley   Input Parameter:
196af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
206af0ca60SMatthew G. Knepley 
216af0ca60SMatthew G. Knepley   Level: beginner
226af0ca60SMatthew G. Knepley 
23db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
246af0ca60SMatthew G. Knepley @*/
259371c9d4SSatish Balay PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) {
266af0ca60SMatthew G. Knepley   PetscFunctionBegin;
276af0ca60SMatthew G. Knepley   if (!*ce) PetscFunctionReturn(0);
286af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific((*ce), PETSC_OBJECT_CLASSID, 1);
296af0ca60SMatthew G. Knepley   if (--((PetscObject)(*ce))->refct > 0) {
306af0ca60SMatthew G. Knepley     *ce = NULL;
316af0ca60SMatthew G. Knepley     PetscFunctionReturn(0);
326af0ca60SMatthew G. Knepley   }
339566063dSJacob Faibussowitsch   PetscCall(PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs));
349566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*ce)->dofs, (*ce)->errors));
359566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(ce));
366af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
376af0ca60SMatthew G. Knepley }
386af0ca60SMatthew G. Knepley 
396af0ca60SMatthew G. Knepley /*@
406af0ca60SMatthew G. Knepley   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
416af0ca60SMatthew G. Knepley 
42*f6dfbefdSBarry Smith   Collective on ce
436af0ca60SMatthew G. Knepley 
446af0ca60SMatthew G. Knepley   Input Parameters:
456af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
466af0ca60SMatthew G. Knepley 
476af0ca60SMatthew G. Knepley   Level: beginner
486af0ca60SMatthew G. Knepley 
49db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
506af0ca60SMatthew G. Knepley @*/
519371c9d4SSatish Balay PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) {
526af0ca60SMatthew G. Knepley   PetscFunctionBegin;
53d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)ce), "", "Convergence Estimator Options", "PetscConvEst");
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL));
58d0609cedSBarry Smith   PetscOptionsEnd();
596af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
606af0ca60SMatthew G. Knepley }
616af0ca60SMatthew G. Knepley 
626af0ca60SMatthew G. Knepley /*@
636af0ca60SMatthew G. Knepley   PetscConvEstView - Views a PetscConvEst object
646af0ca60SMatthew G. Knepley 
65*f6dfbefdSBarry Smith   Collective on ce
666af0ca60SMatthew G. Knepley 
676af0ca60SMatthew G. Knepley   Input Parameters:
686af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
696af0ca60SMatthew G. Knepley - viewer - The PetscViewer object
706af0ca60SMatthew G. Knepley 
716af0ca60SMatthew G. Knepley   Level: beginner
726af0ca60SMatthew G. Knepley 
73db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
746af0ca60SMatthew G. Knepley @*/
759371c9d4SSatish Balay PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) {
766af0ca60SMatthew G. Knepley   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ce, viewer));
7863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %" PetscInt_FMT " levels\n", ce->Nr + 1));
796af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
806af0ca60SMatthew G. Knepley }
816af0ca60SMatthew G. Knepley 
826af0ca60SMatthew G. Knepley /*@
836af0ca60SMatthew G. Knepley   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
846af0ca60SMatthew G. Knepley 
856af0ca60SMatthew G. Knepley   Not collective
866af0ca60SMatthew G. Knepley 
876af0ca60SMatthew G. Knepley   Input Parameter:
886af0ca60SMatthew G. Knepley . ce     - The PetscConvEst object
896af0ca60SMatthew G. Knepley 
906af0ca60SMatthew G. Knepley   Output Parameter:
91900f6b5bSMatthew G. Knepley . solver - The solver
926af0ca60SMatthew G. Knepley 
936af0ca60SMatthew G. Knepley   Level: intermediate
946af0ca60SMatthew G. Knepley 
95db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
966af0ca60SMatthew G. Knepley @*/
979371c9d4SSatish Balay PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver) {
986af0ca60SMatthew G. Knepley   PetscFunctionBegin;
996af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
100900f6b5bSMatthew G. Knepley   PetscValidPointer(solver, 2);
101900f6b5bSMatthew G. Knepley   *solver = ce->solver;
1026af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1036af0ca60SMatthew G. Knepley }
1046af0ca60SMatthew G. Knepley 
1056af0ca60SMatthew G. Knepley /*@
1066af0ca60SMatthew G. Knepley   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
1076af0ca60SMatthew G. Knepley 
1086af0ca60SMatthew G. Knepley   Not collective
1096af0ca60SMatthew G. Knepley 
1106af0ca60SMatthew G. Knepley   Input Parameters:
1116af0ca60SMatthew G. Knepley + ce     - The PetscConvEst object
112900f6b5bSMatthew G. Knepley - solver - The solver
1136af0ca60SMatthew G. Knepley 
1146af0ca60SMatthew G. Knepley   Level: intermediate
1156af0ca60SMatthew G. Knepley 
1166af0ca60SMatthew G. Knepley   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
1176af0ca60SMatthew G. Knepley 
118db781477SPatrick Sanan .seealso: `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1196af0ca60SMatthew G. Knepley @*/
1209371c9d4SSatish Balay PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver) {
1216af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1226af0ca60SMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
123900f6b5bSMatthew G. Knepley   PetscValidHeader(solver, 2);
124900f6b5bSMatthew G. Knepley   ce->solver = solver;
125dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, setsolver, solver);
1266af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1276af0ca60SMatthew G. Knepley }
1286af0ca60SMatthew G. Knepley 
1296af0ca60SMatthew G. Knepley /*@
1306af0ca60SMatthew G. Knepley   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
1316af0ca60SMatthew G. Knepley 
132*f6dfbefdSBarry Smith   Collective on ce
1336af0ca60SMatthew G. Knepley 
1346af0ca60SMatthew G. Knepley   Input Parameters:
1356af0ca60SMatthew G. Knepley . ce - The PetscConvEst object
1366af0ca60SMatthew G. Knepley 
1376af0ca60SMatthew G. Knepley   Level: beginner
1386af0ca60SMatthew G. Knepley 
139db781477SPatrick Sanan .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`
1406af0ca60SMatthew G. Knepley @*/
1419371c9d4SSatish Balay PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) {
142083401c6SMatthew G. Knepley   PetscInt Nf, f, Nds, s;
1436af0ca60SMatthew G. Knepley 
1446af0ca60SMatthew G. Knepley   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(ce->idm, &Nf));
146900f6b5bSMatthew G. Knepley   ce->Nf = PetscMax(Nf, 1);
1479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2((ce->Nr + 1) * ce->Nf, &ce->dofs, (ce->Nr + 1) * ce->Nf, &ce->errors));
1489566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs));
149900f6b5bSMatthew G. Knepley   for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
1509566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(ce->idm, &Nds));
151083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
152083401c6SMatthew G. Knepley     PetscDS         ds;
153083401c6SMatthew G. Knepley     DMLabel         label;
154083401c6SMatthew G. Knepley     IS              fieldIS;
155083401c6SMatthew G. Knepley     const PetscInt *fields;
156083401c6SMatthew G. Knepley     PetscInt        dsNf;
157083401c6SMatthew G. Knepley 
1589566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds));
1599566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &dsNf));
1609566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields));
161083401c6SMatthew G. Knepley     for (f = 0; f < dsNf; ++f) {
162083401c6SMatthew G. Knepley       const PetscInt field = fields[f];
1639566063dSJacob Faibussowitsch       PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]));
164083401c6SMatthew G. Knepley     }
1659566063dSJacob Faibussowitsch     if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields));
166083401c6SMatthew G. Knepley   }
167ad540459SPierre 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);
1686af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1696af0ca60SMatthew G. Knepley }
1706af0ca60SMatthew G. Knepley 
1719371c9d4SSatish Balay PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u) {
1726af0ca60SMatthew G. Knepley   PetscFunctionBegin;
173900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
174900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
175900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
176dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, initguess, r, dm, u);
177900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
178900f6b5bSMatthew G. Knepley }
179900f6b5bSMatthew G. Knepley 
1809371c9d4SSatish Balay PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) {
181900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
182900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
183900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
184900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
185900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
186dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, computeerror, r, dm, u, errors);
187900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
188900f6b5bSMatthew G. Knepley }
189900f6b5bSMatthew G. Knepley 
190f2cacb80SMatthew G. Knepley /*@
191f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
192f2cacb80SMatthew G. Knepley 
193*f6dfbefdSBarry Smith   Collective on ce
194f2cacb80SMatthew G. Knepley 
195d8d19677SJose E. Roman   Input Parameters:
196*f6dfbefdSBarry Smith + ce - The `PetscConvEst` object
197f2cacb80SMatthew G. Knepley - r  - The refinement level
198f2cacb80SMatthew G. Knepley 
199f2cacb80SMatthew G. Knepley   Options database keys:
200ee300463SSatish Balay . -convest_monitor - Activate the monitor
201f2cacb80SMatthew G. Knepley 
202f2cacb80SMatthew G. Knepley   Level: intermediate
203f2cacb80SMatthew G. Knepley 
204*f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
205f2cacb80SMatthew G. Knepley @*/
2069371c9d4SSatish Balay PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r) {
207900f6b5bSMatthew G. Knepley   MPI_Comm comm;
208900f6b5bSMatthew G. Knepley   PetscInt f;
209900f6b5bSMatthew G. Knepley 
210900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
211900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2127809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r * ce->Nf];
213900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r * ce->Nf];
214900f6b5bSMatthew G. Knepley 
2159566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)ce, &comm));
2169566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "N: "));
2179566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
2187809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2199566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
22063a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f]));
2217809adefSMatthew G. Knepley     }
2229566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2239566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "  "));
2249566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Error: "));
2259566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
226900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2279566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2289566063dSJacob Faibussowitsch       if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11"));
2299566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "%g", (double)errors[f]));
230900f6b5bSMatthew G. Knepley     }
2319566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2329566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n"));
233900f6b5bSMatthew G. Knepley   }
234900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
235900f6b5bSMatthew G. Knepley }
236900f6b5bSMatthew G. Knepley 
2379371c9d4SSatish Balay static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver) {
238900f6b5bSMatthew G. Knepley   PetscClassId id;
239900f6b5bSMatthew G. Knepley 
240900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
24208401ef6SPierre Jolivet   PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
2439566063dSJacob Faibussowitsch   PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm));
244900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
245900f6b5bSMatthew G. Knepley }
246900f6b5bSMatthew G. Knepley 
2479371c9d4SSatish Balay static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) {
248900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2499566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u));
250900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
251900f6b5bSMatthew G. Knepley }
252900f6b5bSMatthew G. Knepley 
2539371c9d4SSatish Balay static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) {
254900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors));
256900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
257900f6b5bSMatthew G. Knepley }
258900f6b5bSMatthew G. Knepley 
2599371c9d4SSatish Balay static PetscErrorCode PetscConvEstSetJacobianNullspace_Private(PetscConvEst ce, SNES snes) {
260478db826SMatthew G. Knepley   DM       dm;
261478db826SMatthew G. Knepley   PetscInt f;
262478db826SMatthew G. Knepley 
263478db826SMatthew G. Knepley   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
265478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
266478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
267478db826SMatthew G. Knepley 
2689566063dSJacob Faibussowitsch     PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr));
269478db826SMatthew G. Knepley     if (nspconstr) {
270478db826SMatthew G. Knepley       MatNullSpace nullsp;
271478db826SMatthew G. Knepley       Mat          J;
272478db826SMatthew G. Knepley 
2739566063dSJacob Faibussowitsch       PetscCall((*nspconstr)(dm, f, f, &nullsp));
2749566063dSJacob Faibussowitsch       PetscCall(SNESSetUp(snes));
2759566063dSJacob Faibussowitsch       PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
2769566063dSJacob Faibussowitsch       PetscCall(MatSetNullSpace(J, nullsp));
2779566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullsp));
278478db826SMatthew G. Knepley       break;
279478db826SMatthew G. Knepley     }
280478db826SMatthew G. Knepley   }
281478db826SMatthew G. Knepley   PetscFunctionReturn(0);
282478db826SMatthew G. Knepley }
283478db826SMatthew G. Knepley 
2849371c9d4SSatish Balay static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[]) {
285900f6b5bSMatthew G. Knepley   SNES        snes = (SNES)ce->solver;
286900f6b5bSMatthew G. Knepley   DM         *dm;
287900f6b5bSMatthew G. Knepley   PetscObject disc;
288900f6b5bSMatthew G. Knepley   PetscReal  *x, *y, slope, intercept;
2897809adefSMatthew G. Knepley   PetscInt    Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
290900f6b5bSMatthew G. Knepley   void       *ctx;
291900f6b5bSMatthew G. Knepley 
292900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
29308401ef6SPierre Jolivet   PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r);
2949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
2959566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
2969566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
2979566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
2989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr + 1), &dm));
2996af0ca60SMatthew G. Knepley   /* Loop over meshes */
300900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3016af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
302900f6b5bSMatthew G. Knepley     Vec u;
303e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG)
304adfa7136SMatthew G. Knepley     PetscLogStage stage;
305e5ed2c37SJose E. Roman #endif
306adfa7136SMatthew G. Knepley     char        stageName[PETSC_MAX_PATH_LEN];
307900f6b5bSMatthew G. Knepley     const char *dmname, *uname;
308adfa7136SMatthew G. Knepley 
30963a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r));
310608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
3119566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
3129566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
313608e5a7aSMatthew G. Knepley #endif
3149566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3156af0ca60SMatthew G. Knepley     if (r > 0) {
316b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
3179566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r]));
3189566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r - 1]));
319b2df8587SMatthew G. Knepley       } else {
320b2df8587SMatthew G. Knepley         DM cdm, rcdm;
321b2df8587SMatthew G. Knepley 
3229566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r - 1], &dm[r]));
3239566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r - 1], dm[r]));
3249566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm));
3259566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r], &rcdm));
3269566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
327b2df8587SMatthew G. Knepley       }
3289566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
3299566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname));
3309566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname));
331478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3328cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
333900f6b5bSMatthew G. Knepley 
3349566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr));
3359566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr));
3368f892730SMatthew G. Knepley       }
3376af0ca60SMatthew G. Knepley     }
3389566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
3396af0ca60SMatthew G. Knepley     /* Create solution */
3409566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
3419566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
3429566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
3439566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)u, uname));
3446af0ca60SMatthew G. Knepley     /* Setup solver */
3459566063dSJacob Faibussowitsch     PetscCall(SNESReset(snes));
3469566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm[r]));
3479566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx));
3489566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
349478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
3509566063dSJacob Faibussowitsch     PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
3516af0ca60SMatthew G. Knepley     /* Create initial guess */
3529566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
3539566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
3549566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
3559566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf]));
3569566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
357adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
358a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
359a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
360a56d3bf6SMatthew G. Knepley 
361a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
3629566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
3639566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
3649566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
3659566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
3669566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f]));
3679566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f]));
368adfa7136SMatthew G. Knepley     }
3696af0ca60SMatthew G. Knepley     /* Monitor */
3709566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
37173269098SMatthew G. Knepley     if (!r) {
37273269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
37373269098SMatthew G. Knepley       KSP ksp;
37473269098SMatthew G. Knepley       PC  pc;
37573269098SMatthew G. Knepley 
3769566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
3779566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
3789566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
37973269098SMatthew G. Knepley     }
3806af0ca60SMatthew G. Knepley     /* Cleanup */
3819566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
3829566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
3836af0ca60SMatthew G. Knepley   }
38448a46eb9SPierre Jolivet   for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r]));
3856af0ca60SMatthew G. Knepley   /* Fit convergence rate */
3869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y));
38746079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
3886af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
3897809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]);
39046079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]);
3916af0ca60SMatthew G. Knepley     }
3929566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept));
3936af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
39446079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
39546079b62SMatthew G. Knepley   }
3969566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
3979566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
3982cae373cSMatthew G. Knepley   /* Restore solver */
3999566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
40073269098SMatthew G. Knepley   {
40173269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
40273269098SMatthew G. Knepley     KSP ksp;
40373269098SMatthew G. Knepley     PC  pc;
40473269098SMatthew G. Knepley 
4059566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4069566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4079566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
4089566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
40973269098SMatthew G. Knepley   }
4109566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, ce->idm));
4119566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx));
4129566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4139566063dSJacob Faibussowitsch   PetscCall(PetscConvEstSetJacobianNullspace_Private(ce, snes));
414900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
415900f6b5bSMatthew G. Knepley }
416900f6b5bSMatthew G. Knepley 
417900f6b5bSMatthew G. Knepley /*@
418900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
419900f6b5bSMatthew G. Knepley 
420900f6b5bSMatthew G. Knepley   Not collective
421900f6b5bSMatthew G. Knepley 
422900f6b5bSMatthew G. Knepley   Input Parameter:
423*f6dfbefdSBarry Smith . ce   - The `PetscConvEst` object
424900f6b5bSMatthew G. Knepley 
425900f6b5bSMatthew G. Knepley   Output Parameter:
426900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
427900f6b5bSMatthew G. Knepley 
428*f6dfbefdSBarry Smith   Options database keys:
429*f6dfbefdSBarry Smith + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate
430*f6dfbefdSBarry Smith - -ts_convergence_estimate - Execute convergence estimation inside `TSSolve()` and print out the rate
431*f6dfbefdSBarry Smith 
432*f6dfbefdSBarry Smith   Notes:
433*f6dfbefdSBarry Smith   The convergence rate alpha is defined by
434900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
435900f6b5bSMatthew G. Knepley   where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
436900f6b5bSMatthew G. Knepley   spatial resolution and \Delta t for the temporal resolution.
437900f6b5bSMatthew G. Knepley 
438900f6b5bSMatthew G. Knepley   We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
439900f6b5bSMatthew G. Knepley   based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
440900f6b5bSMatthew G. Knepley 
441900f6b5bSMatthew G. Knepley   Level: intermediate
442900f6b5bSMatthew G. Knepley 
443db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
444900f6b5bSMatthew G. Knepley @*/
4459371c9d4SSatish Balay PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) {
446900f6b5bSMatthew G. Knepley   PetscInt f;
447900f6b5bSMatthew G. Knepley 
448900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
4499566063dSJacob Faibussowitsch   if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event));
450900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
451dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, getconvrate, alpha);
4526af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4536af0ca60SMatthew G. Knepley }
4546af0ca60SMatthew G. Knepley 
4556af0ca60SMatthew G. Knepley /*@
4566af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4576af0ca60SMatthew G. Knepley 
458*f6dfbefdSBarry Smith    Collective on snes
4596af0ca60SMatthew G. Knepley 
4606af0ca60SMatthew G. Knepley    Parameter:
461*f6dfbefdSBarry Smith +  snes - iterative context obtained from `SNESCreate()`
46246079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4636af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4646af0ca60SMatthew G. Knepley 
4656af0ca60SMatthew G. Knepley    Options Database Keys:
4666af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4676af0ca60SMatthew G. Knepley 
4686af0ca60SMatthew G. Knepley    Level: developer
4696af0ca60SMatthew G. Knepley 
470*f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetRate()`
4716af0ca60SMatthew G. Knepley @*/
4729371c9d4SSatish Balay PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) {
4736af0ca60SMatthew G. Knepley   PetscBool isAscii;
4746af0ca60SMatthew G. Knepley 
4756af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
4776af0ca60SMatthew G. Knepley   if (isAscii) {
478900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
479a56d3bf6SMatthew G. Knepley 
4809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel));
4819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: "));
4829566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "["));
483a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
4849566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
4859566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f]));
48646079b62SMatthew G. Knepley     }
4879566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]"));
4889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
4899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel));
4906af0ca60SMatthew G. Knepley   }
4916af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4926af0ca60SMatthew G. Knepley }
493900f6b5bSMatthew G. Knepley 
494900f6b5bSMatthew G. Knepley /*@
495*f6dfbefdSBarry Smith   PetscConvEstCreate - Create a `PetscConvEst` object
496900f6b5bSMatthew G. Knepley 
497900f6b5bSMatthew G. Knepley   Collective
498900f6b5bSMatthew G. Knepley 
499900f6b5bSMatthew G. Knepley   Input Parameter:
500*f6dfbefdSBarry Smith . comm - The communicator for the `PetscConvEst` object
501900f6b5bSMatthew G. Knepley 
502900f6b5bSMatthew G. Knepley   Output Parameter:
503*f6dfbefdSBarry Smith . ce   - The `PetscConvEst` object
504900f6b5bSMatthew G. Knepley 
505900f6b5bSMatthew G. Knepley   Level: beginner
506900f6b5bSMatthew G. Knepley 
507*f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()`
508900f6b5bSMatthew G. Knepley @*/
5099371c9d4SSatish Balay PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) {
510900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
511900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
5129566063dSJacob Faibussowitsch   PetscCall(PetscSysInitializePackage());
5139566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView));
514900f6b5bSMatthew G. Knepley   (*ce)->monitor           = PETSC_FALSE;
5152e61be88SMatthew G. Knepley   (*ce)->r                 = 2.0;
516900f6b5bSMatthew G. Knepley   (*ce)->Nr                = 4;
517900f6b5bSMatthew G. Knepley   (*ce)->event             = -1;
518900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
519900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
520900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
521900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
522900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
523900f6b5bSMatthew G. Knepley }
524