xref: /petsc/src/snes/utils/convest.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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;
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 
17*c3339decSBarry Smith   Collective
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 @*/
26d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
27d71ae5a4SJacob Faibussowitsch {
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 
44*c3339decSBarry Smith   Collective
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 @*/
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();
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 
68*c3339decSBarry Smith   Collective
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 @*/
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));
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 @*/
101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
102d71ae5a4SJacob Faibussowitsch {
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 @*/
125d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
126d71ae5a4SJacob Faibussowitsch {
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;
131dbbe0bcdSBarry 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 
138*c3339decSBarry Smith   Collective
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 @*/
147d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
148d71ae5a4SJacob Faibussowitsch {
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   }
174ad540459SPierre 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);
1756af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
1766af0ca60SMatthew G. Knepley }
1776af0ca60SMatthew G. Knepley 
178d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
179d71ae5a4SJacob Faibussowitsch {
1806af0ca60SMatthew G. Knepley   PetscFunctionBegin;
181900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
182900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
183900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
184dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, initguess, r, dm, u);
185900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
186900f6b5bSMatthew G. Knepley }
187900f6b5bSMatthew G. Knepley 
188d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
189d71ae5a4SJacob Faibussowitsch {
190900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
191900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
192900f6b5bSMatthew G. Knepley   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
193900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
194900f6b5bSMatthew G. Knepley   PetscValidRealPointer(errors, 5);
195dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, computeerror, r, dm, u, errors);
196900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
197900f6b5bSMatthew G. Knepley }
198900f6b5bSMatthew G. Knepley 
199f2cacb80SMatthew G. Knepley /*@
200f2cacb80SMatthew G. Knepley   PetscConvEstMonitorDefault - Monitors the convergence estimation loop
201f2cacb80SMatthew G. Knepley 
202*c3339decSBarry Smith   Collective
203f2cacb80SMatthew G. Knepley 
204d8d19677SJose E. Roman   Input Parameters:
205f6dfbefdSBarry Smith + ce - The `PetscConvEst` object
206f2cacb80SMatthew G. Knepley - r  - The refinement level
207f2cacb80SMatthew G. Knepley 
208f2cacb80SMatthew G. Knepley   Options database keys:
209ee300463SSatish Balay . -convest_monitor - Activate the monitor
210f2cacb80SMatthew G. Knepley 
211f2cacb80SMatthew G. Knepley   Level: intermediate
212f2cacb80SMatthew G. Knepley 
213f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
214f2cacb80SMatthew G. Knepley @*/
215d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
216d71ae5a4SJacob Faibussowitsch {
217900f6b5bSMatthew G. Knepley   MPI_Comm comm;
218900f6b5bSMatthew G. Knepley   PetscInt f;
219900f6b5bSMatthew G. Knepley 
220900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
221900f6b5bSMatthew G. Knepley   if (ce->monitor) {
2227809adefSMatthew G. Knepley     PetscInt  *dofs   = &ce->dofs[r * ce->Nf];
223900f6b5bSMatthew G. Knepley     PetscReal *errors = &ce->errors[r * ce->Nf];
224900f6b5bSMatthew G. Knepley 
2259566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)ce, &comm));
2269566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "N: "));
2279566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
2287809adefSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2299566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
23063a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f]));
2317809adefSMatthew G. Knepley     }
2329566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "  "));
2349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Error: "));
2359566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "["));
236900f6b5bSMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
2379566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscPrintf(comm, ", "));
2389566063dSJacob Faibussowitsch       if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11"));
2399566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "%g", (double)errors[f]));
240900f6b5bSMatthew G. Knepley     }
2419566063dSJacob Faibussowitsch     if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]"));
2429566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n"));
243900f6b5bSMatthew G. Knepley   }
244900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
245900f6b5bSMatthew G. Knepley }
246900f6b5bSMatthew G. Knepley 
247d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
248d71ae5a4SJacob Faibussowitsch {
249900f6b5bSMatthew G. Knepley   PetscClassId id;
250900f6b5bSMatthew G. Knepley 
251900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
25308401ef6SPierre Jolivet   PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
2549566063dSJacob Faibussowitsch   PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm));
255900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
256900f6b5bSMatthew G. Knepley }
257900f6b5bSMatthew G. Knepley 
258d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
259d71ae5a4SJacob Faibussowitsch {
260900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u));
262900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
263900f6b5bSMatthew G. Knepley }
264900f6b5bSMatthew G. Knepley 
265d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
266d71ae5a4SJacob Faibussowitsch {
267900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors));
269900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
270900f6b5bSMatthew G. Knepley }
271900f6b5bSMatthew G. Knepley 
2728c401167SJed Brown static PetscErrorCode PetscConvEstSetJacobianNullSpace_Private(PetscConvEst ce, SNES snes)
273d71ae5a4SJacob Faibussowitsch {
274478db826SMatthew G. Knepley   DM       dm;
275478db826SMatthew G. Knepley   PetscInt f;
276478db826SMatthew G. Knepley 
277478db826SMatthew G. Knepley   PetscFunctionBegin;
2789566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
279478db826SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
280478db826SMatthew G. Knepley     PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
281478db826SMatthew G. Knepley 
2829566063dSJacob Faibussowitsch     PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr));
283478db826SMatthew G. Knepley     if (nspconstr) {
284478db826SMatthew G. Knepley       MatNullSpace nullsp;
285478db826SMatthew G. Knepley       Mat          J;
286478db826SMatthew G. Knepley 
2879566063dSJacob Faibussowitsch       PetscCall((*nspconstr)(dm, f, f, &nullsp));
2889566063dSJacob Faibussowitsch       PetscCall(SNESSetUp(snes));
2899566063dSJacob Faibussowitsch       PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
2909566063dSJacob Faibussowitsch       PetscCall(MatSetNullSpace(J, nullsp));
2919566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullsp));
292478db826SMatthew G. Knepley       break;
293478db826SMatthew G. Knepley     }
294478db826SMatthew G. Knepley   }
295478db826SMatthew G. Knepley   PetscFunctionReturn(0);
296478db826SMatthew G. Knepley }
297478db826SMatthew G. Knepley 
298d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
299d71ae5a4SJacob Faibussowitsch {
300900f6b5bSMatthew G. Knepley   SNES        snes = (SNES)ce->solver;
301900f6b5bSMatthew G. Knepley   DM         *dm;
302900f6b5bSMatthew G. Knepley   PetscObject disc;
303900f6b5bSMatthew G. Knepley   PetscReal  *x, *y, slope, intercept;
3047809adefSMatthew G. Knepley   PetscInt    Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
305900f6b5bSMatthew G. Knepley   void       *ctx;
306900f6b5bSMatthew G. Knepley 
307900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
30808401ef6SPierre Jolivet   PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r);
3099566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
3109566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
3119566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
3129566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
3139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr + 1), &dm));
3146af0ca60SMatthew G. Knepley   /* Loop over meshes */
315900f6b5bSMatthew G. Knepley   dm[0] = ce->idm;
3166af0ca60SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
317900f6b5bSMatthew G. Knepley     Vec u;
318e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG)
319adfa7136SMatthew G. Knepley     PetscLogStage stage;
320e5ed2c37SJose E. Roman #endif
321adfa7136SMatthew G. Knepley     char        stageName[PETSC_MAX_PATH_LEN];
322900f6b5bSMatthew G. Knepley     const char *dmname, *uname;
323adfa7136SMatthew G. Knepley 
32463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r));
325608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
3269566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
3279566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
328608e5a7aSMatthew G. Knepley #endif
3299566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3306af0ca60SMatthew G. Knepley     if (r > 0) {
331b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
3329566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r]));
3339566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r - 1]));
334b2df8587SMatthew G. Knepley       } else {
335b2df8587SMatthew G. Knepley         DM cdm, rcdm;
336b2df8587SMatthew G. Knepley 
3379566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r - 1], &dm[r]));
3389566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r - 1], dm[r]));
3399566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm));
3409566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r], &rcdm));
3419566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
342b2df8587SMatthew G. Knepley       }
3439566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
3449566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname));
3459566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname));
346478db826SMatthew G. Knepley       for (f = 0; f < ce->Nf; ++f) {
3478cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
348900f6b5bSMatthew G. Knepley 
3499566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr));
3509566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr));
3518f892730SMatthew G. Knepley       }
3526af0ca60SMatthew G. Knepley     }
3539566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
3546af0ca60SMatthew G. Knepley     /* Create solution */
3559566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
3569566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
3579566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
3589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)u, uname));
3596af0ca60SMatthew G. Knepley     /* Setup solver */
3609566063dSJacob Faibussowitsch     PetscCall(SNESReset(snes));
3619566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm[r]));
3629566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx));
3639566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
364478db826SMatthew G. Knepley     /* Set nullspace for Jacobian */
3658c401167SJed Brown     PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes));
3666af0ca60SMatthew G. Knepley     /* Create initial guess */
3679566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
3689566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
3699566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
3709566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf]));
3719566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
372adfa7136SMatthew G. Knepley     for (f = 0; f < ce->Nf; ++f) {
373a56d3bf6SMatthew G. Knepley       PetscSection s, fs;
374a56d3bf6SMatthew G. Knepley       PetscInt     lsize;
375a56d3bf6SMatthew G. Knepley 
376a56d3bf6SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
3779566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
3789566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
3799566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
3809566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
3819566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f]));
3829566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f]));
383adfa7136SMatthew G. Knepley     }
3846af0ca60SMatthew G. Knepley     /* Monitor */
3859566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
38673269098SMatthew G. Knepley     if (!r) {
38773269098SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
38873269098SMatthew G. Knepley       KSP ksp;
38973269098SMatthew G. Knepley       PC  pc;
39073269098SMatthew G. Knepley 
3919566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
3929566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
3939566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
39473269098SMatthew G. Knepley     }
3956af0ca60SMatthew G. Knepley     /* Cleanup */
3969566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
3979566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
3986af0ca60SMatthew G. Knepley   }
39948a46eb9SPierre Jolivet   for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r]));
4006af0ca60SMatthew G. Knepley   /* Fit convergence rate */
4019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y));
40246079b62SMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) {
4036af0ca60SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
4047809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]);
40546079b62SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]);
4066af0ca60SMatthew G. Knepley     }
4079566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept));
4086af0ca60SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
40946079b62SMatthew G. Knepley     alpha[f] = -slope * dim;
41046079b62SMatthew G. Knepley   }
4119566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
4132cae373cSMatthew G. Knepley   /* Restore solver */
4149566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
41573269098SMatthew G. Knepley   {
41673269098SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
41773269098SMatthew G. Knepley     KSP ksp;
41873269098SMatthew G. Knepley     PC  pc;
41973269098SMatthew G. Knepley 
4209566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4219566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4229566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
4239566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
42473269098SMatthew G. Knepley   }
4259566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, ce->idm));
4269566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx));
4279566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4288c401167SJed Brown   PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes));
429900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
430900f6b5bSMatthew G. Knepley }
431900f6b5bSMatthew G. Knepley 
432900f6b5bSMatthew G. Knepley /*@
433900f6b5bSMatthew G. Knepley   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
434900f6b5bSMatthew G. Knepley 
435900f6b5bSMatthew G. Knepley   Not collective
436900f6b5bSMatthew G. Knepley 
437900f6b5bSMatthew G. Knepley   Input Parameter:
438f6dfbefdSBarry Smith . ce   - The `PetscConvEst` object
439900f6b5bSMatthew G. Knepley 
440900f6b5bSMatthew G. Knepley   Output Parameter:
441900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field
442900f6b5bSMatthew G. Knepley 
443f6dfbefdSBarry Smith   Options database keys:
444f6dfbefdSBarry Smith + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate
445f6dfbefdSBarry Smith - -ts_convergence_estimate - Execute convergence estimation inside `TSSolve()` and print out the rate
446f6dfbefdSBarry Smith 
447f6dfbefdSBarry Smith   Notes:
448f6dfbefdSBarry Smith   The convergence rate alpha is defined by
449900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha
450900f6b5bSMatthew G. Knepley   where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
451900f6b5bSMatthew G. Knepley   spatial resolution and \Delta t for the temporal resolution.
452900f6b5bSMatthew G. Knepley 
453900f6b5bSMatthew G. Knepley   We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
454900f6b5bSMatthew G. Knepley   based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
455900f6b5bSMatthew G. Knepley 
456900f6b5bSMatthew G. Knepley   Level: intermediate
457900f6b5bSMatthew G. Knepley 
458db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`
459900f6b5bSMatthew G. Knepley @*/
460d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
461d71ae5a4SJacob Faibussowitsch {
462900f6b5bSMatthew G. Knepley   PetscInt f;
463900f6b5bSMatthew G. Knepley 
464900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event));
466900f6b5bSMatthew G. Knepley   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
467dbbe0bcdSBarry Smith   PetscUseTypeMethod(ce, getconvrate, alpha);
4686af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
4696af0ca60SMatthew G. Knepley }
4706af0ca60SMatthew G. Knepley 
4716af0ca60SMatthew G. Knepley /*@
4726af0ca60SMatthew G. Knepley   PetscConvEstRateView - Displays the convergence rate to a viewer
4736af0ca60SMatthew G. Knepley 
474*c3339decSBarry Smith    Collective
4756af0ca60SMatthew G. Knepley 
4766af0ca60SMatthew G. Knepley    Parameter:
477f6dfbefdSBarry Smith +  snes - iterative context obtained from `SNESCreate()`
47846079b62SMatthew G. Knepley .  alpha - the convergence rate for each field
4796af0ca60SMatthew G. Knepley -  viewer - the viewer to display the reason
4806af0ca60SMatthew G. Knepley 
4816af0ca60SMatthew G. Knepley    Options Database Keys:
4826af0ca60SMatthew G. Knepley .  -snes_convergence_estimate - print the convergence rate
4836af0ca60SMatthew G. Knepley 
4846af0ca60SMatthew G. Knepley    Level: developer
4856af0ca60SMatthew G. Knepley 
486f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetRate()`
4876af0ca60SMatthew G. Knepley @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
489d71ae5a4SJacob Faibussowitsch {
4906af0ca60SMatthew G. Knepley   PetscBool isAscii;
4916af0ca60SMatthew G. Knepley 
4926af0ca60SMatthew G. Knepley   PetscFunctionBegin;
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
4946af0ca60SMatthew G. Knepley   if (isAscii) {
495900f6b5bSMatthew G. Knepley     PetscInt Nf = ce->Nf, f;
496a56d3bf6SMatthew G. Knepley 
4979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel));
4989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: "));
4999566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "["));
500a56d3bf6SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5019566063dSJacob Faibussowitsch       if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
5029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f]));
50346079b62SMatthew G. Knepley     }
5049566063dSJacob Faibussowitsch     if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]"));
5059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
5069566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel));
5076af0ca60SMatthew G. Knepley   }
5086af0ca60SMatthew G. Knepley   PetscFunctionReturn(0);
5096af0ca60SMatthew G. Knepley }
510900f6b5bSMatthew G. Knepley 
511900f6b5bSMatthew G. Knepley /*@
512f6dfbefdSBarry Smith   PetscConvEstCreate - Create a `PetscConvEst` object
513900f6b5bSMatthew G. Knepley 
514900f6b5bSMatthew G. Knepley   Collective
515900f6b5bSMatthew G. Knepley 
516900f6b5bSMatthew G. Knepley   Input Parameter:
517f6dfbefdSBarry Smith . comm - The communicator for the `PetscConvEst` object
518900f6b5bSMatthew G. Knepley 
519900f6b5bSMatthew G. Knepley   Output Parameter:
520f6dfbefdSBarry Smith . ce   - The `PetscConvEst` object
521900f6b5bSMatthew G. Knepley 
522900f6b5bSMatthew G. Knepley   Level: beginner
523900f6b5bSMatthew G. Knepley 
524f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()`
525900f6b5bSMatthew G. Knepley @*/
526d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
527d71ae5a4SJacob Faibussowitsch {
528900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
529900f6b5bSMatthew G. Knepley   PetscValidPointer(ce, 2);
5309566063dSJacob Faibussowitsch   PetscCall(PetscSysInitializePackage());
5319566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView));
532900f6b5bSMatthew G. Knepley   (*ce)->monitor           = PETSC_FALSE;
5332e61be88SMatthew G. Knepley   (*ce)->r                 = 2.0;
534900f6b5bSMatthew G. Knepley   (*ce)->Nr                = 4;
535900f6b5bSMatthew G. Knepley   (*ce)->event             = -1;
536900f6b5bSMatthew G. Knepley   (*ce)->ops->setsolver    = PetscConvEstSetSNES_Private;
537900f6b5bSMatthew G. Knepley   (*ce)->ops->initguess    = PetscConvEstInitGuessSNES_Private;
538900f6b5bSMatthew G. Knepley   (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
539900f6b5bSMatthew G. Knepley   (*ce)->ops->getconvrate  = PetscConvEstGetConvRateSNES_Private;
540900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
541900f6b5bSMatthew G. Knepley }
542