xref: /petsc/src/snes/utils/dmplexsnes.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>   /*I "petscsnes.h"   I*/
324cdb843SMatthew G. Knepley #include <petscds.h>
4af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
5af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
6552f7358SJed Brown 
7d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
8d2b2dc1eSMatthew G. Knepley   #include <petscdmceed.h>
9d2b2dc1eSMatthew G. Knepley   #include <petscdmplexceed.h>
10d2b2dc1eSMatthew G. Knepley #endif
11d2b2dc1eSMatthew G. Knepley 
12d71ae5a4SJacob Faibussowitsch static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
13d71ae5a4SJacob Faibussowitsch {
14649ef022SMatthew Knepley   p[0] = u[uOff[1]];
15649ef022SMatthew Knepley }
16649ef022SMatthew Knepley 
17649ef022SMatthew Knepley /*
1820f4b53cSBarry Smith   SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero.
1920f4b53cSBarry Smith   This is normally used only to evaluate convergence rates for the pressure accurately.
20649ef022SMatthew Knepley 
21c3339decSBarry Smith   Collective
22649ef022SMatthew Knepley 
23649ef022SMatthew Knepley   Input Parameters:
24420bcc1bSBarry Smith + snes      - The `SNES`
25649ef022SMatthew Knepley . pfield    - The field number for pressure
26649ef022SMatthew Knepley . nullspace - The pressure nullspace
27649ef022SMatthew Knepley . u         - The solution vector
28649ef022SMatthew Knepley - ctx       - An optional user context
29649ef022SMatthew Knepley 
30649ef022SMatthew Knepley   Output Parameter:
31649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
32649ef022SMatthew Knepley 
332fe279fdSBarry Smith   Level: developer
342fe279fdSBarry Smith 
35420bcc1bSBarry Smith   Note:
36649ef022SMatthew Knepley   If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly.
37649ef022SMatthew Knepley 
38420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESConvergedCorrectPressure()`
39649ef022SMatthew Knepley */
40d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
41d71ae5a4SJacob Faibussowitsch {
42649ef022SMatthew Knepley   DM          dm;
43649ef022SMatthew Knepley   PetscDS     ds;
44649ef022SMatthew Knepley   const Vec  *nullvecs;
45649ef022SMatthew Knepley   PetscScalar pintd, *intc, *intn;
46649ef022SMatthew Knepley   MPI_Comm    comm;
47649ef022SMatthew Knepley   PetscInt    Nf, Nv;
48649ef022SMatthew Knepley 
49649ef022SMatthew Knepley   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
519566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
5228b400f6SJacob Faibussowitsch   PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
5328b400f6SJacob Faibussowitsch   PetscCheck(nullspace, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
549566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
559566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, pfield, pressure_Private));
569566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs));
5763a3b9bcSJacob Faibussowitsch   PetscCheck(Nv == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %" PetscInt_FMT, Nv);
589566063dSJacob Faibussowitsch   PetscCall(VecDot(nullvecs[0], u, &pintd));
5908401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double)PetscRealPart(pintd));
609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &intc, Nf, &intn));
629566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u, -intc[pfield] / intn[pfield], nullvecs[0]));
65649ef022SMatthew Knepley #if defined(PETSC_USE_DEBUG)
669566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
6708401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(intc[pfield]) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double)PetscRealPart(intc[pfield]));
68649ef022SMatthew Knepley #endif
699566063dSJacob Faibussowitsch   PetscCall(PetscFree2(intc, intn));
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71649ef022SMatthew Knepley }
72649ef022SMatthew Knepley 
73649ef022SMatthew Knepley /*@C
74ceaaa498SBarry Smith   SNESConvergedCorrectPressure - The regular `SNES` convergence test that, up on convergence, adds a vector in the nullspace
75ceaaa498SBarry Smith   to make the continuum integral of the pressure field equal to zero.
76649ef022SMatthew Knepley 
77c3339decSBarry Smith   Logically Collective
78649ef022SMatthew Knepley 
79649ef022SMatthew Knepley   Input Parameters:
80f6dfbefdSBarry Smith + snes  - the `SNES` context
81649ef022SMatthew Knepley . it    - the iteration (0 indicates before any Newton steps)
82649ef022SMatthew Knepley . xnorm - 2-norm of current iterate
83e4094ef1SJacob Faibussowitsch . gnorm - 2-norm of current step
84e4094ef1SJacob Faibussowitsch . f     - 2-norm of function at current iterate
85649ef022SMatthew Knepley - ctx   - Optional user context
86649ef022SMatthew Knepley 
87649ef022SMatthew Knepley   Output Parameter:
88f6dfbefdSBarry Smith . reason - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FNORM_NAN`
89649ef022SMatthew Knepley 
9020f4b53cSBarry Smith   Options Database Key:
9120f4b53cSBarry Smith . -snes_convergence_test correct_pressure - see `SNESSetFromOptions()`
9220f4b53cSBarry Smith 
9320f4b53cSBarry Smith   Level: advanced
9420f4b53cSBarry Smith 
95649ef022SMatthew Knepley   Notes:
96f6dfbefdSBarry Smith   In order to use this convergence test, you must set up several PETSc structures. First fields must be added to the `DM`, and a `PetscDS`
97f6dfbefdSBarry Smith   must be created with discretizations of those fields. We currently assume that the pressure field has index 1.
98f6dfbefdSBarry Smith   The pressure field must have a nullspace, likely created using the `DMSetNullSpaceConstructor()` interface.
99f6dfbefdSBarry Smith   Last we must be able to integrate the pressure over the domain, so the `DM` attached to the SNES `must` be a `DMPLEX` at this time.
100f362779dSJed Brown 
101ceaaa498SBarry Smith   Developer Note:
102ceaaa498SBarry Smith   This is a total misuse of the `SNES` convergence test handling system. It should be removed. Perhaps a `SNESSetPostSolve()` could
103ceaaa498SBarry Smith   be constructed to handle this process.
104ceaaa498SBarry Smith 
105420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DM`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()`
106649ef022SMatthew Knepley @*/
107d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
108d71ae5a4SJacob Faibussowitsch {
109649ef022SMatthew Knepley   PetscBool monitorIntegral = PETSC_FALSE;
110649ef022SMatthew Knepley 
111649ef022SMatthew Knepley   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx));
113649ef022SMatthew Knepley   if (monitorIntegral) {
114649ef022SMatthew Knepley     Mat          J;
115649ef022SMatthew Knepley     Vec          u;
116649ef022SMatthew Knepley     MatNullSpace nullspace;
117649ef022SMatthew Knepley     const Vec   *nullvecs;
118649ef022SMatthew Knepley     PetscScalar  pintd;
119649ef022SMatthew Knepley 
1209566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1219566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1229566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1239566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs));
1249566063dSJacob Faibussowitsch     PetscCall(VecDot(nullvecs[0], u, &pintd));
1259566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double)PetscRealPart(pintd)));
126649ef022SMatthew Knepley   }
127649ef022SMatthew Knepley   if (*reason > 0) {
128649ef022SMatthew Knepley     Mat          J;
129649ef022SMatthew Knepley     Vec          u;
130649ef022SMatthew Knepley     MatNullSpace nullspace;
131649ef022SMatthew Knepley     PetscInt     pfield = 1;
132649ef022SMatthew Knepley 
1339566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1349566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1359566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1369566063dSJacob Faibussowitsch     PetscCall(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx));
137649ef022SMatthew Knepley   }
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139649ef022SMatthew Knepley }
140649ef022SMatthew Knepley 
141d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
142d71ae5a4SJacob Faibussowitsch {
1436da023fcSToby Isaac   PetscBool isPlex;
1446da023fcSToby Isaac 
1456da023fcSToby Isaac   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
1476da023fcSToby Isaac   if (isPlex) {
1486da023fcSToby Isaac     *plex = dm;
1499566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
150f7148743SMatthew G. Knepley   } else {
1519566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
152f7148743SMatthew G. Knepley     if (!*plex) {
1539566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
1549566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
155cbf8eb3cSStefano Zampini     } else {
156cbf8eb3cSStefano Zampini       PetscCall(PetscObjectReference((PetscObject)*plex));
157cbf8eb3cSStefano Zampini     }
1586da023fcSToby Isaac     if (copy) {
1599566063dSJacob Faibussowitsch       PetscCall(DMCopyDMSNES(dm, *plex));
1609566063dSJacob Faibussowitsch       PetscCall(DMCopyAuxiliaryVec(dm, *plex));
1616da023fcSToby Isaac     }
1626da023fcSToby Isaac   }
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1646da023fcSToby Isaac }
1656da023fcSToby Isaac 
1664267b1a3SMatthew G. Knepley /*@C
167cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
168cc0c4584SMatthew G. Knepley 
169c3339decSBarry Smith   Collective
170cc0c4584SMatthew G. Knepley 
171cc0c4584SMatthew G. Knepley   Input Parameters:
172420bcc1bSBarry Smith + snes   - the `SNES` context, must have an attached `DM`
173cc0c4584SMatthew G. Knepley . its    - iteration number
174cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
175f6dfbefdSBarry Smith - vf     - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII`
176cc0c4584SMatthew G. Knepley 
1772fe279fdSBarry Smith   Level: intermediate
1782fe279fdSBarry Smith 
179f6dfbefdSBarry Smith   Note:
180cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
181cc0c4584SMatthew G. Knepley 
182420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()`
183cc0c4584SMatthew G. Knepley @*/
184d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
185d71ae5a4SJacob Faibussowitsch {
186d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
187cc0c4584SMatthew G. Knepley   Vec                res;
188cc0c4584SMatthew G. Knepley   DM                 dm;
189cc0c4584SMatthew G. Knepley   PetscSection       s;
190cc0c4584SMatthew G. Knepley   const PetscScalar *r;
191cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
192cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
193cc0c4584SMatthew G. Knepley 
194cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
1954d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
1969566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
1979566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1989566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
1999566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
2009566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
2019566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms));
2029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(res, &r));
203cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
204cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
205cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
206cc0c4584SMatthew G. Knepley 
2079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
2089566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
209cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d]));
210cc0c4584SMatthew G. Knepley     }
211cc0c4584SMatthew G. Knepley   }
2129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(res, &r));
213*e91c04dfSPierre Jolivet   PetscCallMPI(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
2149566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
2159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
21663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm));
217cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
2189566063dSJacob Faibussowitsch     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
2199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f])));
220cc0c4584SMatthew G. Knepley   }
2219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
2229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
2239566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
2249566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lnorms, norms));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226cc0c4584SMatthew G. Knepley }
22724cdb843SMatthew G. Knepley 
2286493148fSStefano Zampini /********************* SNES callbacks **************************/
2296493148fSStefano Zampini 
2306493148fSStefano Zampini /*@
2316493148fSStefano Zampini   DMPlexSNESComputeObjectiveFEM - Sums the local objectives from the local input X using pointwise functions specified by the user
2326493148fSStefano Zampini 
2336493148fSStefano Zampini   Input Parameters:
2346493148fSStefano Zampini + dm   - The mesh
2356493148fSStefano Zampini . X    - Local solution
2366493148fSStefano Zampini - user - The user context
2376493148fSStefano Zampini 
2386493148fSStefano Zampini   Output Parameter:
2396493148fSStefano Zampini . obj - Local objective value
2406493148fSStefano Zampini 
2416493148fSStefano Zampini   Level: developer
2426493148fSStefano Zampini 
2436493148fSStefano Zampini .seealso: `DM`, `DMPlexSNESComputeResidualFEM()`
2446493148fSStefano Zampini @*/
2456493148fSStefano Zampini PetscErrorCode DMPlexSNESComputeObjectiveFEM(DM dm, Vec X, PetscReal *obj, void *user)
2466493148fSStefano Zampini {
2476493148fSStefano Zampini   PetscInt     Nf, cellHeight, cStart, cEnd;
2486493148fSStefano Zampini   PetscScalar *cintegral;
2496493148fSStefano Zampini 
2506493148fSStefano Zampini   PetscFunctionBegin;
2516493148fSStefano Zampini   PetscCall(DMGetNumFields(dm, &Nf));
2526493148fSStefano Zampini   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2536493148fSStefano Zampini   PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
2546493148fSStefano Zampini   PetscCall(PetscCalloc1((cEnd - cStart) * Nf, &cintegral));
2556493148fSStefano Zampini   PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2566493148fSStefano Zampini   PetscCall(DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user));
2576493148fSStefano Zampini   /* Sum up values */
2586493148fSStefano Zampini   *obj = 0;
2596493148fSStefano Zampini   for (PetscInt c = cStart; c < cEnd; ++c)
2606493148fSStefano Zampini     for (PetscInt f = 0; f < Nf; ++f) *obj += PetscRealPart(cintegral[(c - cStart) * Nf + f]);
2616493148fSStefano Zampini   PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2626493148fSStefano Zampini   PetscCall(PetscFree(cintegral));
2636493148fSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2646493148fSStefano Zampini }
26524cdb843SMatthew G. Knepley 
26624cdb843SMatthew G. Knepley /*@
267420bcc1bSBarry Smith   DMPlexSNESComputeResidualFEM - Sums the local residual into vector `F` from the local input `X` using pointwise functions specified by the user
26824cdb843SMatthew G. Knepley 
26924cdb843SMatthew G. Knepley   Input Parameters:
27024cdb843SMatthew G. Knepley + dm   - The mesh
27124cdb843SMatthew G. Knepley . X    - Local solution
27224cdb843SMatthew G. Knepley - user - The user context
27324cdb843SMatthew G. Knepley 
27424cdb843SMatthew G. Knepley   Output Parameter:
27524cdb843SMatthew G. Knepley . F - Local output vector
27624cdb843SMatthew G. Knepley 
2772fe279fdSBarry Smith   Level: developer
2782fe279fdSBarry Smith 
279f6dfbefdSBarry Smith   Note:
280420bcc1bSBarry Smith   The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional.
2818db23a53SJed Brown 
2826493148fSStefano Zampini .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMSNESComputeJacobianAction()`
28324cdb843SMatthew G. Knepley @*/
284d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
285d71ae5a4SJacob Faibussowitsch {
2866da023fcSToby Isaac   DM       plex;
287083401c6SMatthew G. Knepley   IS       allcellIS;
2886528b96dSMatthew G. Knepley   PetscInt Nds, s;
28924cdb843SMatthew G. Knepley 
29024cdb843SMatthew G. Knepley   PetscFunctionBegin;
2919566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
2929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2939566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
2946528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
2956528b96dSMatthew G. Knepley     PetscDS      ds;
2966528b96dSMatthew G. Knepley     IS           cellIS;
29706ad1575SMatthew G. Knepley     PetscFormKey key;
2986528b96dSMatthew G. Knepley 
29907218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
3006528b96dSMatthew G. Knepley     key.value = 0;
3016528b96dSMatthew G. Knepley     key.field = 0;
30206ad1575SMatthew G. Knepley     key.part  = 0;
3036528b96dSMatthew G. Knepley     if (!key.label) {
3049566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
3056528b96dSMatthew G. Knepley       cellIS = allcellIS;
3066528b96dSMatthew G. Knepley     } else {
3076528b96dSMatthew G. Knepley       IS pointIS;
3086528b96dSMatthew G. Knepley 
3096528b96dSMatthew G. Knepley       key.value = 1;
3109566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
3119566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
3129566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
3136528b96dSMatthew G. Knepley     }
3149566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
3159566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
3166528b96dSMatthew G. Knepley   }
3179566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
3189566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3206528b96dSMatthew G. Knepley }
3216528b96dSMatthew G. Knepley 
322bdd6f66aSToby Isaac /*@
323420bcc1bSBarry Smith   DMPlexSNESComputeResidualDS - Sums the local residual into vector `F` from the local input `X` using all pointwise functions with unique keys in the `PetscDS`
324d2b2dc1eSMatthew G. Knepley 
325d2b2dc1eSMatthew G. Knepley   Input Parameters:
326d2b2dc1eSMatthew G. Knepley + dm   - The mesh
327d2b2dc1eSMatthew G. Knepley . X    - Local solution
328d2b2dc1eSMatthew G. Knepley - user - The user context
329d2b2dc1eSMatthew G. Knepley 
330d2b2dc1eSMatthew G. Knepley   Output Parameter:
331d2b2dc1eSMatthew G. Knepley . F - Local output vector
332d2b2dc1eSMatthew G. Knepley 
333d2b2dc1eSMatthew G. Knepley   Level: developer
334d2b2dc1eSMatthew G. Knepley 
335d2b2dc1eSMatthew G. Knepley   Note:
336420bcc1bSBarry Smith   The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional.
337d2b2dc1eSMatthew G. Knepley 
338420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()`
339d2b2dc1eSMatthew G. Knepley @*/
340d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualDS(DM dm, Vec X, Vec F, void *user)
341d2b2dc1eSMatthew G. Knepley {
342d2b2dc1eSMatthew G. Knepley   DM       plex;
343d2b2dc1eSMatthew G. Knepley   IS       allcellIS;
344d2b2dc1eSMatthew G. Knepley   PetscInt Nds, s;
345d2b2dc1eSMatthew G. Knepley 
346d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
347d2b2dc1eSMatthew G. Knepley   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
348d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
349d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetNumDS(dm, &Nds));
350d2b2dc1eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
351d2b2dc1eSMatthew G. Knepley     PetscDS ds;
352d2b2dc1eSMatthew G. Knepley     DMLabel label;
353d2b2dc1eSMatthew G. Knepley     IS      cellIS;
354d2b2dc1eSMatthew G. Knepley 
355d2b2dc1eSMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL));
356d2b2dc1eSMatthew G. Knepley     {
357d2b2dc1eSMatthew G. Knepley       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
358d2b2dc1eSMatthew G. Knepley       PetscWeakForm     wf;
359d2b2dc1eSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
360d2b2dc1eSMatthew G. Knepley       PetscFormKey     *reskeys;
361d2b2dc1eSMatthew G. Knepley 
362d2b2dc1eSMatthew G. Knepley       /* Get unique residual keys */
363d2b2dc1eSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
364d2b2dc1eSMatthew G. Knepley         PetscInt Nkm;
365d2b2dc1eSMatthew G. Knepley         PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm));
366d2b2dc1eSMatthew G. Knepley         Nk += Nkm;
367d2b2dc1eSMatthew G. Knepley       }
368d2b2dc1eSMatthew G. Knepley       PetscCall(PetscMalloc1(Nk, &reskeys));
369d2b2dc1eSMatthew G. Knepley       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys));
370d2b2dc1eSMatthew G. Knepley       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
371d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFormKeySort(Nk, reskeys));
372d2b2dc1eSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
373d2b2dc1eSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
374d2b2dc1eSMatthew G. Knepley           ++k;
375d2b2dc1eSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
376d2b2dc1eSMatthew G. Knepley         }
377d2b2dc1eSMatthew G. Knepley       }
378d2b2dc1eSMatthew G. Knepley       Nk = k;
379d2b2dc1eSMatthew G. Knepley 
380d2b2dc1eSMatthew G. Knepley       PetscCall(PetscDSGetWeakForm(ds, &wf));
381d2b2dc1eSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
382d2b2dc1eSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
383d2b2dc1eSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
384d2b2dc1eSMatthew G. Knepley 
385d2b2dc1eSMatthew G. Knepley         if (!label) {
386d2b2dc1eSMatthew G. Knepley           PetscCall(PetscObjectReference((PetscObject)allcellIS));
387d2b2dc1eSMatthew G. Knepley           cellIS = allcellIS;
388d2b2dc1eSMatthew G. Knepley         } else {
389d2b2dc1eSMatthew G. Knepley           IS pointIS;
390d2b2dc1eSMatthew G. Knepley 
391d2b2dc1eSMatthew G. Knepley           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
392d2b2dc1eSMatthew G. Knepley           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
393d2b2dc1eSMatthew G. Knepley           PetscCall(ISDestroy(&pointIS));
394d2b2dc1eSMatthew G. Knepley         }
395d2b2dc1eSMatthew G. Knepley         PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
396d2b2dc1eSMatthew G. Knepley         PetscCall(ISDestroy(&cellIS));
397d2b2dc1eSMatthew G. Knepley       }
398d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFree(reskeys));
399d2b2dc1eSMatthew G. Knepley     }
400d2b2dc1eSMatthew G. Knepley   }
401d2b2dc1eSMatthew G. Knepley   PetscCall(ISDestroy(&allcellIS));
402d2b2dc1eSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
403d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
404d2b2dc1eSMatthew G. Knepley }
405d2b2dc1eSMatthew G. Knepley 
406d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
407d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualCEED(DM dm, Vec locX, Vec locF, void *user)
408d2b2dc1eSMatthew G. Knepley {
409d2b2dc1eSMatthew G. Knepley   Ceed       ceed;
410d2b2dc1eSMatthew G. Knepley   DMCeed     sd = dm->dmceed;
411d2b2dc1eSMatthew G. Knepley   CeedVector clocX, clocF;
412d2b2dc1eSMatthew G. Knepley 
413d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
414d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
415d2b2dc1eSMatthew G. Knepley   PetscCheck(sd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DM has no CEED data. Call DMCeedCreate() before computing the residual.");
416d2b2dc1eSMatthew G. Knepley   PetscCall(DMCeedComputeGeometry(dm, sd));
417d2b2dc1eSMatthew G. Knepley 
418d2b2dc1eSMatthew G. Knepley   PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX));
419d2b2dc1eSMatthew G. Knepley   PetscCall(VecGetCeedVector(locF, ceed, &clocF));
420d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE));
421d2b2dc1eSMatthew G. Knepley   PetscCall(VecRestoreCeedVectorRead(locX, &clocX));
422d2b2dc1eSMatthew G. Knepley   PetscCall(VecRestoreCeedVector(locF, &clocF));
423d2b2dc1eSMatthew G. Knepley 
424d2b2dc1eSMatthew G. Knepley   {
425d2b2dc1eSMatthew G. Knepley     DM_Plex *mesh = (DM_Plex *)dm->data;
426d2b2dc1eSMatthew G. Knepley 
427d2b2dc1eSMatthew G. Knepley     if (mesh->printFEM) {
428d2b2dc1eSMatthew G. Knepley       PetscSection section;
429d2b2dc1eSMatthew G. Knepley       Vec          locFbc;
430d2b2dc1eSMatthew G. Knepley       PetscInt     pStart, pEnd, p, maxDof;
431d2b2dc1eSMatthew G. Knepley       PetscScalar *zeroes;
432d2b2dc1eSMatthew G. Knepley 
433d2b2dc1eSMatthew G. Knepley       PetscCall(DMGetLocalSection(dm, &section));
434d2b2dc1eSMatthew G. Knepley       PetscCall(VecDuplicate(locF, &locFbc));
435d2b2dc1eSMatthew G. Knepley       PetscCall(VecCopy(locF, locFbc));
436d2b2dc1eSMatthew G. Knepley       PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
437d2b2dc1eSMatthew G. Knepley       PetscCall(PetscSectionGetMaxDof(section, &maxDof));
438d2b2dc1eSMatthew G. Knepley       PetscCall(PetscCalloc1(maxDof, &zeroes));
439d2b2dc1eSMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
440d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFree(zeroes));
441d2b2dc1eSMatthew G. Knepley       PetscCall(DMPrintLocalVec(dm, "Residual", mesh->printTol, locFbc));
442d2b2dc1eSMatthew G. Knepley       PetscCall(VecDestroy(&locFbc));
443d2b2dc1eSMatthew G. Knepley     }
444d2b2dc1eSMatthew G. Knepley   }
445d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
446d2b2dc1eSMatthew G. Knepley }
447d2b2dc1eSMatthew G. Knepley #endif
448d2b2dc1eSMatthew G. Knepley 
449d2b2dc1eSMatthew G. Knepley /*@
450420bcc1bSBarry Smith   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input `X`
451bdd6f66aSToby Isaac 
452bdd6f66aSToby Isaac   Input Parameters:
453bdd6f66aSToby Isaac + dm   - The mesh
454bdd6f66aSToby Isaac - user - The user context
455bdd6f66aSToby Isaac 
456bdd6f66aSToby Isaac   Output Parameter:
457bdd6f66aSToby Isaac . X - Local solution
458bdd6f66aSToby Isaac 
459bdd6f66aSToby Isaac   Level: developer
460bdd6f66aSToby Isaac 
461420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()`
462bdd6f66aSToby Isaac @*/
463d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
464d71ae5a4SJacob Faibussowitsch {
465bdd6f66aSToby Isaac   DM plex;
466bdd6f66aSToby Isaac 
467bdd6f66aSToby Isaac   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
4699566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL));
4709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472bdd6f66aSToby Isaac }
473bdd6f66aSToby Isaac 
4747a73cf09SMatthew G. Knepley /*@
475c118d280SBarry Smith   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(`X`) on `Y`
4767a73cf09SMatthew G. Knepley 
4777a73cf09SMatthew G. Knepley   Input Parameters:
478f6dfbefdSBarry Smith + dm   - The `DM`
4797a73cf09SMatthew G. Knepley . X    - Local solution vector
4807a73cf09SMatthew G. Knepley . Y    - Local input vector
4817a73cf09SMatthew G. Knepley - user - The user context
4827a73cf09SMatthew G. Knepley 
4837a73cf09SMatthew G. Knepley   Output Parameter:
48469d47153SPierre Jolivet . F - local output vector
4857a73cf09SMatthew G. Knepley 
4867a73cf09SMatthew G. Knepley   Level: developer
4877a73cf09SMatthew G. Knepley 
488c118d280SBarry Smith   Note:
489f6dfbefdSBarry Smith   Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly.
4908e3a2eefSMatthew G. Knepley 
491c118d280SBarry Smith   This only works with `DMPLEX`
492c118d280SBarry Smith 
493c118d280SBarry Smith   Developer Note:
494c118d280SBarry Smith   This should be called `DMPlexSNESComputeJacobianAction()`
495c118d280SBarry Smith 
496a94f484eSPierre Jolivet .seealso: [](ch_snes), `DM`, `DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()`
4977a73cf09SMatthew G. Knepley @*/
498d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
499d71ae5a4SJacob Faibussowitsch {
5008e3a2eefSMatthew G. Knepley   DM       plex;
5018e3a2eefSMatthew G. Knepley   IS       allcellIS;
5028e3a2eefSMatthew G. Knepley   PetscInt Nds, s;
503a925c78cSMatthew G. Knepley 
504a925c78cSMatthew G. Knepley   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
5069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
5079566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
5088e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
5098e3a2eefSMatthew G. Knepley     PetscDS ds;
5108e3a2eefSMatthew G. Knepley     DMLabel label;
5118e3a2eefSMatthew G. Knepley     IS      cellIS;
5127a73cf09SMatthew G. Knepley 
51307218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL));
5148e3a2eefSMatthew G. Knepley     {
51506ad1575SMatthew G. Knepley       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
5168e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
5178e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
51806ad1575SMatthew G. Knepley       PetscFormKey     *jackeys;
5198e3a2eefSMatthew G. Knepley 
5208e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
5218e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
5228e3a2eefSMatthew G. Knepley         PetscInt Nkm;
5239566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm));
5248e3a2eefSMatthew G. Knepley         Nk += Nkm;
5258e3a2eefSMatthew G. Knepley       }
5269566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &jackeys));
52748a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys));
52863a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
5299566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, jackeys));
5308e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
5318e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
5328e3a2eefSMatthew G. Knepley           ++k;
5338e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
5348e3a2eefSMatthew G. Knepley         }
5358e3a2eefSMatthew G. Knepley       }
5368e3a2eefSMatthew G. Knepley       Nk = k;
5378e3a2eefSMatthew G. Knepley 
5389566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
5398e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
5408e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
5418e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
5428e3a2eefSMatthew G. Knepley 
5438e3a2eefSMatthew G. Knepley         if (!label) {
5449566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
5458e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
5467a73cf09SMatthew G. Knepley         } else {
5478e3a2eefSMatthew G. Knepley           IS pointIS;
548a925c78cSMatthew G. Knepley 
5499566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
5509566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
5519566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
552a925c78cSMatthew G. Knepley         }
5539566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user));
5549566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
5558e3a2eefSMatthew G. Knepley       }
5569566063dSJacob Faibussowitsch       PetscCall(PetscFree(jackeys));
5578e3a2eefSMatthew G. Knepley     }
5588e3a2eefSMatthew G. Knepley   }
5599566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
5609566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
5613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
562a925c78cSMatthew G. Knepley }
563a925c78cSMatthew G. Knepley 
56424cdb843SMatthew G. Knepley /*@
5652fe279fdSBarry Smith   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix `Jac` at the local solution `X` using pointwise functions specified by the user.
56624cdb843SMatthew G. Knepley 
56724cdb843SMatthew G. Knepley   Input Parameters:
5682fe279fdSBarry Smith + dm   - The `DM`
56924cdb843SMatthew G. Knepley . X    - Local input vector
57024cdb843SMatthew G. Knepley - user - The user context
57124cdb843SMatthew G. Knepley 
5722fe279fdSBarry Smith   Output Parameters:
5732fe279fdSBarry Smith + Jac  - Jacobian matrix
5742fe279fdSBarry Smith - JacP - approximate Jacobian from which the preconditioner will be built, often `Jac`
5752fe279fdSBarry Smith 
5762fe279fdSBarry Smith   Level: developer
57724cdb843SMatthew G. Knepley 
57824cdb843SMatthew G. Knepley   Note:
57924cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
58024cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
58124cdb843SMatthew G. Knepley 
582420bcc1bSBarry Smith .seealso: [](ch_snes), `DMPLEX`, `Mat`
58324cdb843SMatthew G. Knepley @*/
584d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user)
585d71ae5a4SJacob Faibussowitsch {
5866da023fcSToby Isaac   DM        plex;
587083401c6SMatthew G. Knepley   IS        allcellIS;
588f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
5896528b96dSMatthew G. Knepley   PetscInt  Nds, s;
59024cdb843SMatthew G. Knepley 
59124cdb843SMatthew G. Knepley   PetscFunctionBegin;
5929566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
5939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
5949566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
595083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
596083401c6SMatthew G. Knepley     PetscDS      ds;
597083401c6SMatthew G. Knepley     IS           cellIS;
59806ad1575SMatthew G. Knepley     PetscFormKey key;
599083401c6SMatthew G. Knepley 
60007218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
6016528b96dSMatthew G. Knepley     key.value = 0;
6026528b96dSMatthew G. Knepley     key.field = 0;
60306ad1575SMatthew G. Knepley     key.part  = 0;
6046528b96dSMatthew G. Knepley     if (!key.label) {
6059566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
606083401c6SMatthew G. Knepley       cellIS = allcellIS;
607083401c6SMatthew G. Knepley     } else {
608083401c6SMatthew G. Knepley       IS pointIS;
609083401c6SMatthew G. Knepley 
6106528b96dSMatthew G. Knepley       key.value = 1;
6119566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
6129566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
6139566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
614083401c6SMatthew G. Knepley     }
615083401c6SMatthew G. Knepley     if (!s) {
6169566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
6179566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
6189566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
6199566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
620083401c6SMatthew G. Knepley     }
6219566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user));
6229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
623083401c6SMatthew G. Knepley   }
6249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
6259566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62724cdb843SMatthew G. Knepley }
6281878804aSMatthew G. Knepley 
6299371c9d4SSatish Balay struct _DMSNESJacobianMFCtx {
6308e3a2eefSMatthew G. Knepley   DM    dm;
6318e3a2eefSMatthew G. Knepley   Vec   X;
6328e3a2eefSMatthew G. Knepley   void *ctx;
6338e3a2eefSMatthew G. Knepley };
6348e3a2eefSMatthew G. Knepley 
635d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
636d71ae5a4SJacob Faibussowitsch {
6378e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
6388e3a2eefSMatthew G. Knepley 
6398e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
6409566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
6419566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, NULL));
6429566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
6439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->X));
6449566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6468e3a2eefSMatthew G. Knepley }
6478e3a2eefSMatthew G. Knepley 
648d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
649d71ae5a4SJacob Faibussowitsch {
6508e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
6518e3a2eefSMatthew G. Knepley 
6528e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
6539566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
6549566063dSJacob Faibussowitsch   PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx));
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6568e3a2eefSMatthew G. Knepley }
6578e3a2eefSMatthew G. Knepley 
6588e3a2eefSMatthew G. Knepley /*@
659f6dfbefdSBarry Smith   DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free
6608e3a2eefSMatthew G. Knepley 
661c3339decSBarry Smith   Collective
6628e3a2eefSMatthew G. Knepley 
6638e3a2eefSMatthew G. Knepley   Input Parameters:
664f6dfbefdSBarry Smith + dm   - The `DM`
6658e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
6662fe279fdSBarry Smith - user - A user context, or `NULL`
6678e3a2eefSMatthew G. Knepley 
6688e3a2eefSMatthew G. Knepley   Output Parameter:
669f6dfbefdSBarry Smith . J - The `Mat`
6708e3a2eefSMatthew G. Knepley 
6718e3a2eefSMatthew G. Knepley   Level: advanced
6728e3a2eefSMatthew G. Knepley 
673c118d280SBarry Smith   Notes:
6742fe279fdSBarry Smith   Vec `X` is kept in `J`, so updating `X` then updates the evaluation point.
6758e3a2eefSMatthew G. Knepley 
676c118d280SBarry Smith   This only works for `DMPLEX`
677c118d280SBarry Smith 
678420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `SNES`, `DMSNESComputeJacobianAction()`
6798e3a2eefSMatthew G. Knepley @*/
680d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
681d71ae5a4SJacob Faibussowitsch {
6828e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
6838e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
6848e3a2eefSMatthew G. Knepley 
6858e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
6869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
6879566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
6889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
6899566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
6909566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, n, n, N, N));
6919566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
6929566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)X));
6939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &ctx));
6948e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
6958e3a2eefSMatthew G. Knepley   ctx->X   = X;
6968e3a2eefSMatthew G. Knepley   ctx->ctx = user;
6979566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, ctx));
6989566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private));
6999566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private));
7003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7018e3a2eefSMatthew G. Knepley }
7028e3a2eefSMatthew G. Knepley 
703d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
704d71ae5a4SJacob Faibussowitsch {
70538cfc46eSPierre Jolivet   SNES   snes;
70638cfc46eSPierre Jolivet   Mat    pJ;
70738cfc46eSPierre Jolivet   DM     ovldm, origdm;
70838cfc46eSPierre Jolivet   DMSNES sdm;
70938cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM, Vec, void *);
71038cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *);
71138cfc46eSPierre Jolivet   void *bctx, *jctx;
71238cfc46eSPierre Jolivet 
71338cfc46eSPierre Jolivet   PetscFunctionBegin;
7149566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ));
71528b400f6SJacob Faibussowitsch   PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat");
7169566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm));
71728b400f6SJacob Faibussowitsch   PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM");
7189566063dSJacob Faibussowitsch   PetscCall(MatGetDM(pJ, &ovldm));
7199566063dSJacob Faibussowitsch   PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx));
7209566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx));
7219566063dSJacob Faibussowitsch   PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx));
7229566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx));
7239566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes));
72438cfc46eSPierre Jolivet   if (!snes) {
7259566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes));
7269566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, ovldm));
7279566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes));
7289566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)snes));
72938cfc46eSPierre Jolivet   }
7309566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(ovldm, &sdm));
7319566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(X));
732800f99ffSJeremy L Thompson   {
733800f99ffSJeremy L Thompson     void *ctx;
734800f99ffSJeremy L Thompson     PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *);
735800f99ffSJeremy L Thompson     PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx));
736800f99ffSJeremy L Thompson     PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx));
737800f99ffSJeremy L Thompson   }
7389566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(X));
73938cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
74038cfc46eSPierre Jolivet   {
74138cfc46eSPierre Jolivet     Mat locpJ;
74238cfc46eSPierre Jolivet 
7439566063dSJacob Faibussowitsch     PetscCall(MatISGetLocalMat(pJ, &locpJ));
7449566063dSJacob Faibussowitsch     PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN));
74538cfc46eSPierre Jolivet   }
7463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74738cfc46eSPierre Jolivet }
74838cfc46eSPierre Jolivet 
749a925c78cSMatthew G. Knepley /*@
7506493148fSStefano Zampini   DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, objective, residual, and Jacobian.
7519f520fc2SToby Isaac 
7529f520fc2SToby Isaac   Input Parameters:
753f6dfbefdSBarry Smith + dm      - The `DM` object
7546493148fSStefano Zampini . use_obj - Use the objective function callback
7556493148fSStefano Zampini - ctx     - The user context that will be passed to pointwise evaluation routines
7561a244344SSatish Balay 
7571a244344SSatish Balay   Level: developer
758f6dfbefdSBarry Smith 
7596493148fSStefano Zampini .seealso: [](ch_snes),`DMPLEX`, `SNES`, `PetscDSAddBoundary()`, `PetscDSSetObjective()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()`
7609f520fc2SToby Isaac @*/
7616493148fSStefano Zampini PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, PetscBool use_obj, void *ctx)
762d71ae5a4SJacob Faibussowitsch {
763d2b2dc1eSMatthew G. Knepley   PetscBool useCeed;
764d2b2dc1eSMatthew G. Knepley 
7659f520fc2SToby Isaac   PetscFunctionBegin;
766d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
7676493148fSStefano Zampini   PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, ctx));
7686493148fSStefano Zampini   if (use_obj) PetscCall(DMSNESSetObjectiveLocal(dm, DMPlexSNESComputeObjectiveFEM, ctx));
769d2b2dc1eSMatthew G. Knepley   if (useCeed) {
770d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
7716493148fSStefano Zampini     PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualCEED, ctx));
772d2b2dc1eSMatthew G. Knepley #else
773d2b2dc1eSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot use CEED traversals without LibCEED. Rerun configure with --download-ceed");
774d2b2dc1eSMatthew G. Knepley #endif
7756493148fSStefano Zampini   } else PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, ctx));
7766493148fSStefano Zampini   PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, ctx));
7779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex));
7783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7799f520fc2SToby Isaac }
7809f520fc2SToby Isaac 
781cc4c1da9SBarry Smith /*@
782f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
783f017bcb6SMatthew G. Knepley 
784f017bcb6SMatthew G. Knepley   Input Parameters:
785f6dfbefdSBarry Smith + snes - the `SNES` object
786f6dfbefdSBarry Smith . dm   - the `DM`
787f2cacb80SMatthew G. Knepley . t    - the time
788f6dfbefdSBarry Smith . u    - a `DM` vector
789f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
790f017bcb6SMatthew G. Knepley 
7912fe279fdSBarry Smith   Output Parameter:
7922fe279fdSBarry Smith . error - An array which holds the discretization error in each field, or `NULL`
7932fe279fdSBarry Smith 
7942fe279fdSBarry Smith   Level: developer
795f3c94560SMatthew G. Knepley 
796f6dfbefdSBarry Smith   Note:
797f6dfbefdSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
7987f96f943SMatthew G. Knepley 
799420bcc1bSBarry Smith   Developer Note:
800420bcc1bSBarry Smith   How is this related to `PetscConvEst`?
801420bcc1bSBarry Smith 
802420bcc1bSBarry Smith .seealso: [](ch_snes), `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`
803f017bcb6SMatthew G. Knepley @*/
804d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
805d71ae5a4SJacob Faibussowitsch {
806f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
807f017bcb6SMatthew G. Knepley   void     **ectxs;
808f3c94560SMatthew G. Knepley   PetscReal *err;
8097f96f943SMatthew G. Knepley   MPI_Comm   comm;
8107f96f943SMatthew G. Knepley   PetscInt   Nf, f;
8111878804aSMatthew G. Knepley 
8121878804aSMatthew G. Knepley   PetscFunctionBegin;
813f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
814f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
815064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
8164f572ea9SToby Isaac   if (error) PetscAssertPointer(error, 6);
8177f96f943SMatthew G. Knepley 
8189566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
8199566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
8207f96f943SMatthew G. Knepley 
8219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
8229566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
8239566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
8247f96f943SMatthew G. Knepley   {
8257f96f943SMatthew G. Knepley     PetscInt Nds, s;
8267f96f943SMatthew G. Knepley 
8279566063dSJacob Faibussowitsch     PetscCall(DMGetNumDS(dm, &Nds));
828083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
8297f96f943SMatthew G. Knepley       PetscDS         ds;
830083401c6SMatthew G. Knepley       DMLabel         label;
831083401c6SMatthew G. Knepley       IS              fieldIS;
8327f96f943SMatthew G. Knepley       const PetscInt *fields;
8337f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
834083401c6SMatthew G. Knepley 
83507218a29SMatthew G. Knepley       PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
8369566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(ds, &dsNf));
8379566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
838083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
839083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
8409566063dSJacob Faibussowitsch         PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
841083401c6SMatthew G. Knepley       }
8429566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
843083401c6SMatthew G. Knepley     }
844083401c6SMatthew G. Knepley   }
845f017bcb6SMatthew G. Knepley   if (Nf > 1) {
8469566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err));
847f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
848ad540459SPierre Jolivet       for (f = 0; f < Nf; ++f) PetscCheck(err[f] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %" PetscInt_FMT " exceeds tolerance %g", (double)err[f], f, (double)tol);
849b878532fSMatthew G. Knepley     } else if (error) {
850b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
8511878804aSMatthew G. Knepley     } else {
8529566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: ["));
853f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
8549566063dSJacob Faibussowitsch         if (f) PetscCall(PetscPrintf(comm, ", "));
8559566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(comm, "%g", (double)err[f]));
8561878804aSMatthew G. Knepley       }
8579566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "]\n"));
858f017bcb6SMatthew G. Knepley     }
859f017bcb6SMatthew G. Knepley   } else {
8609566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]));
861f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
86208401ef6SPierre Jolivet       PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol);
863b878532fSMatthew G. Knepley     } else if (error) {
864b878532fSMatthew G. Knepley       error[0] = err[0];
865f017bcb6SMatthew G. Knepley     } else {
8669566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]));
867f017bcb6SMatthew G. Knepley     }
868f017bcb6SMatthew G. Knepley   }
8699566063dSJacob Faibussowitsch   PetscCall(PetscFree3(exacts, ectxs, err));
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871f017bcb6SMatthew G. Knepley }
872f017bcb6SMatthew G. Knepley 
873cc4c1da9SBarry Smith /*@
874f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
875f017bcb6SMatthew G. Knepley 
876f017bcb6SMatthew G. Knepley   Input Parameters:
877f6dfbefdSBarry Smith + snes - the `SNES` object
878f6dfbefdSBarry Smith . dm   - the `DM`
879f6dfbefdSBarry Smith . u    - a `DM` vector
880f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
881f017bcb6SMatthew G. Knepley 
882f6dfbefdSBarry Smith   Output Parameter:
8832fe279fdSBarry Smith . residual - The residual norm of the exact solution, or `NULL`
884f3c94560SMatthew G. Knepley 
885f017bcb6SMatthew G. Knepley   Level: developer
886f017bcb6SMatthew G. Knepley 
887420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
888f017bcb6SMatthew G. Knepley @*/
889d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
890d71ae5a4SJacob Faibussowitsch {
891f017bcb6SMatthew G. Knepley   MPI_Comm  comm;
892f017bcb6SMatthew G. Knepley   Vec       r;
893f017bcb6SMatthew G. Knepley   PetscReal res;
894f017bcb6SMatthew G. Knepley 
895f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
896f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
897f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
898f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
8994f572ea9SToby Isaac   if (residual) PetscAssertPointer(residual, 5);
9009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
9019566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
9029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
9039566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, u, r));
9049566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
905f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
90608401ef6SPierre Jolivet     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
907b878532fSMatthew G. Knepley   } else if (residual) {
908b878532fSMatthew G. Knepley     *residual = res;
909f017bcb6SMatthew G. Knepley   } else {
9109566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
91193ec0da9SPierre Jolivet     PetscCall(VecFilter(r, 1.0e-10));
9129566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
9139566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
9149566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
915f017bcb6SMatthew G. Knepley   }
9169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
918f017bcb6SMatthew G. Knepley }
919f017bcb6SMatthew G. Knepley 
920cc4c1da9SBarry Smith /*@
921f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
922f017bcb6SMatthew G. Knepley 
923f017bcb6SMatthew G. Knepley   Input Parameters:
924f6dfbefdSBarry Smith + snes - the `SNES` object
925f6dfbefdSBarry Smith . dm   - the `DM`
926f6dfbefdSBarry Smith . u    - a `DM` vector
927f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
928f017bcb6SMatthew G. Knepley 
929f3c94560SMatthew G. Knepley   Output Parameters:
9302fe279fdSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL`
9312fe279fdSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL`
932f3c94560SMatthew G. Knepley 
933f017bcb6SMatthew G. Knepley   Level: developer
934f017bcb6SMatthew G. Knepley 
935420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
936f017bcb6SMatthew G. Knepley @*/
937d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
938d71ae5a4SJacob Faibussowitsch {
939f017bcb6SMatthew G. Knepley   MPI_Comm     comm;
940f017bcb6SMatthew G. Knepley   PetscDS      ds;
941f017bcb6SMatthew G. Knepley   Mat          J, M;
942f017bcb6SMatthew G. Knepley   MatNullSpace nullspace;
943f3c94560SMatthew G. Knepley   PetscReal    slope, intercept;
9447f96f943SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
945f017bcb6SMatthew G. Knepley 
946f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
947f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
948d5d53b22SStefano Zampini   if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
949d5d53b22SStefano Zampini   if (u) PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
9504f572ea9SToby Isaac   if (isLinear) PetscAssertPointer(isLinear, 5);
9514f572ea9SToby Isaac   if (convRate) PetscAssertPointer(convRate, 6);
9529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
953d5d53b22SStefano Zampini   if (!dm) PetscCall(SNESGetDM(snes, &dm));
954d5d53b22SStefano Zampini   if (u) PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
955d5d53b22SStefano Zampini   else PetscCall(SNESGetSolution(snes, &u));
956f017bcb6SMatthew G. Knepley   /* Create and view matrices */
9579566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
9589566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
9599566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
9609566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
961282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
9629566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
9639566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, M));
9649566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
9659566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
9669566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
9679566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
968282e7bb4SMatthew G. Knepley   } else {
9699566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, J));
970282e7bb4SMatthew G. Knepley   }
9719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
9729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
9739566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
974f017bcb6SMatthew G. Knepley   /* Check nullspace */
9759566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
976f017bcb6SMatthew G. Knepley   if (nullspace) {
9771878804aSMatthew G. Knepley     PetscBool isNull;
9789566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
97928b400f6SJacob Faibussowitsch     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
9801878804aSMatthew G. Knepley   }
981f017bcb6SMatthew G. Knepley   /* Taylor test */
982f017bcb6SMatthew G. Knepley   {
983f017bcb6SMatthew G. Knepley     PetscRandom rand;
984f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
985f3c94560SMatthew G. Knepley     PetscReal   h;
986f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
987f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
988f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
989f017bcb6SMatthew G. Knepley 
990f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
9919566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
9929566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
9939566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
9949566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
9959566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
9969566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
997f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
9989566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
9999566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, u, r));
1000f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1001fbccb6d4SPierre Jolivet     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
10029566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
10039566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
10049566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
1005f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
10069566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
1007f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
10089566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, uhat, rhat));
10099566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
10109566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
1011f017bcb6SMatthew G. Knepley 
101238d69901SMatthew G. Knepley       es[Nv] = errors[Nv] == 0 ? -16.0 : PetscLog10Real(errors[Nv]);
1013f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1014f017bcb6SMatthew G. Knepley     }
10159566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
10169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
10179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
10189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
10199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
1020f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1021f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1022f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
1023f017bcb6SMatthew G. Knepley     }
1024f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
10259566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
10269566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
1027f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1028f017bcb6SMatthew G. Knepley     if (tol >= 0) {
10290b121fc5SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
1030b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1031b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1032b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1033f017bcb6SMatthew G. Knepley     } else {
10349566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
10359566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
1036f017bcb6SMatthew G. Knepley     }
1037f017bcb6SMatthew G. Knepley   }
10389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
10393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10401878804aSMatthew G. Knepley }
10411878804aSMatthew G. Knepley 
1042d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1043d71ae5a4SJacob Faibussowitsch {
1044f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
10459566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL));
10469566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL));
10479566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL));
10483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1049f017bcb6SMatthew G. Knepley }
1050f017bcb6SMatthew G. Knepley 
1051cc4c1da9SBarry Smith /*@
1052bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1053bee9a294SMatthew G. Knepley 
1054bee9a294SMatthew G. Knepley   Input Parameters:
1055f6dfbefdSBarry Smith + snes - the `SNES` object
1056f6dfbefdSBarry Smith - u    - representative `SNES` vector
10577f96f943SMatthew G. Knepley 
10582fe279fdSBarry Smith   Level: developer
10592fe279fdSBarry Smith 
1060f6dfbefdSBarry Smith   Note:
1061420bcc1bSBarry Smith   The user must call `PetscDSSetExactSolution()` before this call
1062bee9a294SMatthew G. Knepley 
1063420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DM`
1064bee9a294SMatthew G. Knepley @*/
1065d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
1066d71ae5a4SJacob Faibussowitsch {
10671878804aSMatthew G. Knepley   DM        dm;
10681878804aSMatthew G. Knepley   Vec       sol;
10691878804aSMatthew G. Knepley   PetscBool check;
10701878804aSMatthew G. Knepley 
10711878804aSMatthew G. Knepley   PetscFunctionBegin;
10729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check));
10733ba16761SJacob Faibussowitsch   if (!check) PetscFunctionReturn(PETSC_SUCCESS);
10749566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
10759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
10769566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, sol));
10779566063dSJacob Faibussowitsch   PetscCall(DMSNESCheck_Internal(snes, dm, sol));
10789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
10793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1080552f7358SJed Brown }
1081