xref: /petsc/src/snes/utils/dmplexsnes.c (revision c3af198d9aa427cfd901dc43379d1993d8858a15)
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:
24649ef022SMatthew Knepley + 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 
35649ef022SMatthew Knepley   Notes:
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 
38db781477SPatrick Sanan .seealso: `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 
105e4094ef1SJacob Faibussowitsch .seealso: `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 
14124cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
14224cdb843SMatthew G. Knepley 
143d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
144d71ae5a4SJacob Faibussowitsch {
1456da023fcSToby Isaac   PetscBool isPlex;
1466da023fcSToby Isaac 
1476da023fcSToby Isaac   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
1496da023fcSToby Isaac   if (isPlex) {
1506da023fcSToby Isaac     *plex = dm;
1519566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
152f7148743SMatthew G. Knepley   } else {
1539566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
154f7148743SMatthew G. Knepley     if (!*plex) {
1559566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
1569566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
1576da023fcSToby Isaac       if (copy) {
1589566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
1599566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
1606da023fcSToby Isaac       }
161f7148743SMatthew G. Knepley     } else {
1629566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*plex));
163f7148743SMatthew G. Knepley     }
1646da023fcSToby Isaac   }
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1666da023fcSToby Isaac }
1676da023fcSToby Isaac 
1684267b1a3SMatthew G. Knepley /*@C
169f6dfbefdSBarry Smith   DMInterpolationCreate - Creates a `DMInterpolationInfo` context
1704267b1a3SMatthew G. Knepley 
171d083f849SBarry Smith   Collective
1724267b1a3SMatthew G. Knepley 
1734267b1a3SMatthew G. Knepley   Input Parameter:
1744267b1a3SMatthew G. Knepley . comm - the communicator
1754267b1a3SMatthew G. Knepley 
1764267b1a3SMatthew G. Knepley   Output Parameter:
1774267b1a3SMatthew G. Knepley . ctx - the context
1784267b1a3SMatthew G. Knepley 
1794267b1a3SMatthew G. Knepley   Level: beginner
1804267b1a3SMatthew G. Knepley 
181f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()`
1824267b1a3SMatthew G. Knepley @*/
183d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
184d71ae5a4SJacob Faibussowitsch {
185552f7358SJed Brown   PetscFunctionBegin;
1864f572ea9SToby Isaac   PetscAssertPointer(ctx, 2);
1879566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
1881aa26658SKarl Rupp 
189552f7358SJed Brown   (*ctx)->comm   = comm;
190552f7358SJed Brown   (*ctx)->dim    = -1;
191552f7358SJed Brown   (*ctx)->nInput = 0;
1920298fd71SBarry Smith   (*ctx)->points = NULL;
1930298fd71SBarry Smith   (*ctx)->cells  = NULL;
194552f7358SJed Brown   (*ctx)->n      = -1;
1950298fd71SBarry Smith   (*ctx)->coords = NULL;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197552f7358SJed Brown }
198552f7358SJed Brown 
1994267b1a3SMatthew G. Knepley /*@C
2004267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
2014267b1a3SMatthew G. Knepley 
20220f4b53cSBarry Smith   Not Collective
2034267b1a3SMatthew G. Knepley 
2044267b1a3SMatthew G. Knepley   Input Parameters:
2054267b1a3SMatthew G. Knepley + ctx - the context
2064267b1a3SMatthew G. Knepley - dim - the spatial dimension
2074267b1a3SMatthew G. Knepley 
2084267b1a3SMatthew G. Knepley   Level: intermediate
2094267b1a3SMatthew G. Knepley 
210f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2114267b1a3SMatthew G. Knepley @*/
212d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
213d71ae5a4SJacob Faibussowitsch {
214552f7358SJed Brown   PetscFunctionBegin;
21563a3b9bcSJacob Faibussowitsch   PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim);
216552f7358SJed Brown   ctx->dim = dim;
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218552f7358SJed Brown }
219552f7358SJed Brown 
2204267b1a3SMatthew G. Knepley /*@C
2214267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2224267b1a3SMatthew G. Knepley 
22320f4b53cSBarry Smith   Not Collective
2244267b1a3SMatthew G. Knepley 
2254267b1a3SMatthew G. Knepley   Input Parameter:
2264267b1a3SMatthew G. Knepley . ctx - the context
2274267b1a3SMatthew G. Knepley 
2284267b1a3SMatthew G. Knepley   Output Parameter:
2294267b1a3SMatthew G. Knepley . dim - the spatial dimension
2304267b1a3SMatthew G. Knepley 
2314267b1a3SMatthew G. Knepley   Level: intermediate
2324267b1a3SMatthew G. Knepley 
233f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2344267b1a3SMatthew G. Knepley @*/
235d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
236d71ae5a4SJacob Faibussowitsch {
237552f7358SJed Brown   PetscFunctionBegin;
2384f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
239552f7358SJed Brown   *dim = ctx->dim;
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241552f7358SJed Brown }
242552f7358SJed Brown 
2434267b1a3SMatthew G. Knepley /*@C
2444267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2454267b1a3SMatthew G. Knepley 
24620f4b53cSBarry Smith   Not Collective
2474267b1a3SMatthew G. Knepley 
2484267b1a3SMatthew G. Knepley   Input Parameters:
2494267b1a3SMatthew G. Knepley + ctx - the context
2504267b1a3SMatthew G. Knepley - dof - the number of fields
2514267b1a3SMatthew G. Knepley 
2524267b1a3SMatthew G. Knepley   Level: intermediate
2534267b1a3SMatthew G. Knepley 
254f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2554267b1a3SMatthew G. Knepley @*/
256d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
257d71ae5a4SJacob Faibussowitsch {
258552f7358SJed Brown   PetscFunctionBegin;
25963a3b9bcSJacob Faibussowitsch   PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof);
260552f7358SJed Brown   ctx->dof = dof;
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262552f7358SJed Brown }
263552f7358SJed Brown 
2644267b1a3SMatthew G. Knepley /*@C
2654267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2664267b1a3SMatthew G. Knepley 
26720f4b53cSBarry Smith   Not Collective
2684267b1a3SMatthew G. Knepley 
2694267b1a3SMatthew G. Knepley   Input Parameter:
2704267b1a3SMatthew G. Knepley . ctx - the context
2714267b1a3SMatthew G. Knepley 
2724267b1a3SMatthew G. Knepley   Output Parameter:
2734267b1a3SMatthew G. Knepley . dof - the number of fields
2744267b1a3SMatthew G. Knepley 
2754267b1a3SMatthew G. Knepley   Level: intermediate
2764267b1a3SMatthew G. Knepley 
277e4094ef1SJacob Faibussowitsch .seealso: `DMInterpolationInfo`, `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2784267b1a3SMatthew G. Knepley @*/
279d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
280d71ae5a4SJacob Faibussowitsch {
281552f7358SJed Brown   PetscFunctionBegin;
2824f572ea9SToby Isaac   PetscAssertPointer(dof, 2);
283552f7358SJed Brown   *dof = ctx->dof;
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
285552f7358SJed Brown }
286552f7358SJed Brown 
2874267b1a3SMatthew G. Knepley /*@C
2884267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2894267b1a3SMatthew G. Knepley 
29020f4b53cSBarry Smith   Not Collective
2914267b1a3SMatthew G. Knepley 
2924267b1a3SMatthew G. Knepley   Input Parameters:
2934267b1a3SMatthew G. Knepley + ctx    - the context
2944267b1a3SMatthew G. Knepley . n      - the number of points
2954267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2964267b1a3SMatthew G. Knepley 
2972fe279fdSBarry Smith   Level: intermediate
2982fe279fdSBarry Smith 
299f6dfbefdSBarry Smith   Note:
300f6dfbefdSBarry Smith   The coordinate information is copied.
3014267b1a3SMatthew G. Knepley 
302f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()`
3034267b1a3SMatthew G. Knepley @*/
304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
305d71ae5a4SJacob Faibussowitsch {
306552f7358SJed Brown   PetscFunctionBegin;
30708401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
30828b400f6SJacob Faibussowitsch   PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
309552f7358SJed Brown   ctx->nInput = n;
3101aa26658SKarl Rupp 
3119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points));
3129566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim));
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
314552f7358SJed Brown }
315552f7358SJed Brown 
3164267b1a3SMatthew G. Knepley /*@C
31752aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
3184267b1a3SMatthew G. Knepley 
319c3339decSBarry Smith   Collective
3204267b1a3SMatthew G. Knepley 
3214267b1a3SMatthew G. Knepley   Input Parameters:
3224267b1a3SMatthew G. Knepley + ctx                 - the context
323f6dfbefdSBarry Smith . dm                  - the `DM` for the function space used for interpolation
324f6dfbefdSBarry Smith . redundantPoints     - If `PETSC_TRUE`, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
325f6dfbefdSBarry Smith - ignoreOutsideDomain - If `PETSC_TRUE`, ignore points outside the domain, otherwise return an error
3264267b1a3SMatthew G. Knepley 
3274267b1a3SMatthew G. Knepley   Level: intermediate
3284267b1a3SMatthew G. Knepley 
329e4094ef1SJacob Faibussowitsch .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
3304267b1a3SMatthew G. Knepley @*/
331d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
332d71ae5a4SJacob Faibussowitsch {
333552f7358SJed Brown   MPI_Comm           comm = ctx->comm;
334552f7358SJed Brown   PetscScalar       *a;
335552f7358SJed Brown   PetscInt           p, q, i;
336552f7358SJed Brown   PetscMPIInt        rank, size;
337552f7358SJed Brown   Vec                pointVec;
3383a93e3b7SToby Isaac   PetscSF            cellSF;
339552f7358SJed Brown   PetscLayout        layout;
340552f7358SJed Brown   PetscReal         *globalPoints;
341cb313848SJed Brown   PetscScalar       *globalPointsScalar;
342552f7358SJed Brown   const PetscInt    *ranges;
343552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3443a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3453a93e3b7SToby Isaac   const PetscInt    *foundPoints;
346552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3473a93e3b7SToby Isaac   PetscInt           n, N, numFound;
348552f7358SJed Brown 
34919436ca2SJed Brown   PetscFunctionBegin;
350064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
3519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
35308401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
35419436ca2SJed Brown   /* Locate points */
35519436ca2SJed Brown   n = ctx->nInput;
356552f7358SJed Brown   if (!redundantPoints) {
3579566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &layout));
3589566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(layout, 1));
3599566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(layout, n));
3609566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(layout));
3619566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(layout, &N));
362552f7358SJed Brown     /* Communicate all points to all processes */
3639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs));
3649566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(layout, &ranges));
365552f7358SJed Brown     for (p = 0; p < size; ++p) {
366552f7358SJed Brown       counts[p] = (ranges[p + 1] - ranges[p]) * ctx->dim;
367552f7358SJed Brown       displs[p] = ranges[p] * ctx->dim;
368552f7358SJed Brown     }
3699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(ctx->points, n * ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm));
370552f7358SJed Brown   } else {
371552f7358SJed Brown     N            = n;
372552f7358SJed Brown     globalPoints = ctx->points;
37338ea73c8SJed Brown     counts = displs = NULL;
37438ea73c8SJed Brown     layout          = NULL;
375552f7358SJed Brown   }
376552f7358SJed Brown #if 0
3779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs));
37819436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
379552f7358SJed Brown #else
380cb313848SJed Brown   #if defined(PETSC_USE_COMPLEX)
3819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar));
38246b3086cSToby Isaac   for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
383cb313848SJed Brown   #else
384cb313848SJed Brown   globalPointsScalar = globalPoints;
385cb313848SJed Brown   #endif
3869566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec));
3879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs));
388ad540459SPierre Jolivet   for (p = 0; p < N; ++p) foundProcs[p] = size;
3893a93e3b7SToby Isaac   cellSF = NULL;
3909566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF));
3919566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells));
392552f7358SJed Brown #endif
3933a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3943a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
395552f7358SJed Brown   }
396552f7358SJed Brown   /* Let the lowest rank process own each point */
3971c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm));
398552f7358SJed Brown   ctx->n = 0;
399552f7358SJed Brown   for (p = 0; p < N; ++p) {
40052aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
4019371c9d4SSatish Balay       PetscCheck(ignoreOutsideDomain, comm, PETSC_ERR_PLIB, "Point %" PetscInt_FMT ": %g %g %g not located in mesh", p, (double)globalPoints[p * ctx->dim + 0], (double)(ctx->dim > 1 ? globalPoints[p * ctx->dim + 1] : 0.0),
4029371c9d4SSatish Balay                  (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0.0));
403f7d195e4SLawrence Mitchell       if (rank == 0) ++ctx->n;
40452aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
405552f7358SJed Brown   }
406552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
4079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->n, &ctx->cells));
4089566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &ctx->coords));
4099566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE));
4109566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(ctx->coords, ctx->dim));
4119566063dSJacob Faibussowitsch   PetscCall(VecSetType(ctx->coords, VECSTANDARD));
4129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->coords, &a));
413552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
414552f7358SJed Brown     if (globalProcs[p] == rank) {
415552f7358SJed Brown       PetscInt d;
416552f7358SJed Brown 
4171aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d];
418f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
419f317b747SMatthew G. Knepley       ++q;
420552f7358SJed Brown     }
421dd400576SPatrick Sanan     if (globalProcs[p] == size && rank == 0) {
42252aa1562SMatthew G. Knepley       PetscInt d;
42352aa1562SMatthew G. Knepley 
42452aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
42552aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
42652aa1562SMatthew G. Knepley       ++q;
42752aa1562SMatthew G. Knepley     }
428552f7358SJed Brown   }
4299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->coords, &a));
430552f7358SJed Brown #if 0
4319566063dSJacob Faibussowitsch   PetscCall(PetscFree3(foundCells,foundProcs,globalProcs));
432552f7358SJed Brown #else
4339566063dSJacob Faibussowitsch   PetscCall(PetscFree2(foundProcs, globalProcs));
4349566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
4359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pointVec));
436552f7358SJed Brown #endif
4379566063dSJacob Faibussowitsch   if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar));
4389566063dSJacob Faibussowitsch   if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs));
4399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
4403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
441552f7358SJed Brown }
442552f7358SJed Brown 
4434267b1a3SMatthew G. Knepley /*@C
444f6dfbefdSBarry Smith   DMInterpolationGetCoordinates - Gets a `Vec` with the coordinates of each interpolation point
4454267b1a3SMatthew G. Knepley 
446c3339decSBarry Smith   Collective
4474267b1a3SMatthew G. Knepley 
4484267b1a3SMatthew G. Knepley   Input Parameter:
4494267b1a3SMatthew G. Knepley . ctx - the context
4504267b1a3SMatthew G. Knepley 
4514267b1a3SMatthew G. Knepley   Output Parameter:
4524267b1a3SMatthew G. Knepley . coordinates - the coordinates of interpolation points
4534267b1a3SMatthew G. Knepley 
45420f4b53cSBarry Smith   Level: intermediate
45520f4b53cSBarry Smith 
456f6dfbefdSBarry Smith   Note:
457f6dfbefdSBarry Smith   The local vector entries correspond to interpolation points lying on this process, according to the associated `DM`.
458f6dfbefdSBarry Smith   This is a borrowed vector that the user should not destroy.
4594267b1a3SMatthew G. Knepley 
460f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4614267b1a3SMatthew G. Knepley @*/
462d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
463d71ae5a4SJacob Faibussowitsch {
464552f7358SJed Brown   PetscFunctionBegin;
4654f572ea9SToby Isaac   PetscAssertPointer(coordinates, 2);
46628b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
467552f7358SJed Brown   *coordinates = ctx->coords;
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
469552f7358SJed Brown }
470552f7358SJed Brown 
4714267b1a3SMatthew G. Knepley /*@C
472f6dfbefdSBarry Smith   DMInterpolationGetVector - Gets a `Vec` which can hold all the interpolated field values
4734267b1a3SMatthew G. Knepley 
474c3339decSBarry Smith   Collective
4754267b1a3SMatthew G. Knepley 
4764267b1a3SMatthew G. Knepley   Input Parameter:
4774267b1a3SMatthew G. Knepley . ctx - the context
4784267b1a3SMatthew G. Knepley 
4794267b1a3SMatthew G. Knepley   Output Parameter:
4804267b1a3SMatthew G. Knepley . v - a vector capable of holding the interpolated field values
4814267b1a3SMatthew G. Knepley 
4822fe279fdSBarry Smith   Level: intermediate
4832fe279fdSBarry Smith 
484f6dfbefdSBarry Smith   Note:
485f6dfbefdSBarry Smith   This vector should be returned using `DMInterpolationRestoreVector()`.
4864267b1a3SMatthew G. Knepley 
487f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4884267b1a3SMatthew G. Knepley @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
490d71ae5a4SJacob Faibussowitsch {
491552f7358SJed Brown   PetscFunctionBegin;
4924f572ea9SToby Isaac   PetscAssertPointer(v, 2);
49328b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4949566063dSJacob Faibussowitsch   PetscCall(VecCreate(ctx->comm, v));
4959566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE));
4969566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*v, ctx->dof));
4979566063dSJacob Faibussowitsch   PetscCall(VecSetType(*v, VECSTANDARD));
4983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
499552f7358SJed Brown }
500552f7358SJed Brown 
5014267b1a3SMatthew G. Knepley /*@C
502f6dfbefdSBarry Smith   DMInterpolationRestoreVector - Returns a `Vec` which can hold all the interpolated field values
5034267b1a3SMatthew G. Knepley 
504c3339decSBarry Smith   Collective
5054267b1a3SMatthew G. Knepley 
5064267b1a3SMatthew G. Knepley   Input Parameters:
5074267b1a3SMatthew G. Knepley + ctx - the context
5084267b1a3SMatthew G. Knepley - v   - a vector capable of holding the interpolated field values
5094267b1a3SMatthew G. Knepley 
5104267b1a3SMatthew G. Knepley   Level: intermediate
5114267b1a3SMatthew G. Knepley 
512f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
5134267b1a3SMatthew G. Knepley @*/
514d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
515d71ae5a4SJacob Faibussowitsch {
516552f7358SJed Brown   PetscFunctionBegin;
5174f572ea9SToby Isaac   PetscAssertPointer(v, 2);
51828b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
5199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(v));
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
521552f7358SJed Brown }
522552f7358SJed Brown 
523ab72731fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v)
524d71ae5a4SJacob Faibussowitsch {
525ab72731fSMatthew G. Knepley   const PetscInt     c   = ctx->cells[p];
526cfe5744fSMatthew G. Knepley   const PetscInt     dof = ctx->dof;
527ab72731fSMatthew G. Knepley   PetscScalar       *x   = NULL;
528cfe5744fSMatthew G. Knepley   const PetscScalar *coords;
529cfe5744fSMatthew G. Knepley   PetscScalar       *a;
530ab72731fSMatthew G. Knepley   PetscReal          v0, J, invJ, detJ, xir[1];
531ab72731fSMatthew G. Knepley   PetscInt           xSize, comp;
532cfe5744fSMatthew G. Knepley 
533cfe5744fSMatthew G. Knepley   PetscFunctionBegin;
5349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5369566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ));
537cfe5744fSMatthew G. Knepley   PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
538cfe5744fSMatthew G. Knepley   xir[0] = invJ * PetscRealPart(coords[p] - v0);
5399566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
540cfe5744fSMatthew G. Knepley   if (2 * dof == xSize) {
541cfe5744fSMatthew G. Knepley     for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) + x[1 * dof + comp] * xir[0];
542cfe5744fSMatthew G. Knepley   } else if (dof == xSize) {
543cfe5744fSMatthew G. Knepley     for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
544cfe5744fSMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input closure size %" PetscInt_FMT " must be either %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 2 * dof, dof);
5459566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
5469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549cfe5744fSMatthew G. Knepley }
550cfe5744fSMatthew G. Knepley 
551ab72731fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v)
552d71ae5a4SJacob Faibussowitsch {
553ab72731fSMatthew G. Knepley   const PetscInt     c = ctx->cells[p];
554ab72731fSMatthew G. Knepley   PetscScalar       *x = NULL;
55556044e6dSMatthew G. Knepley   const PetscScalar *coords;
55656044e6dSMatthew G. Knepley   PetscScalar       *a;
557ab72731fSMatthew G. Knepley   PetscReal         *v0, *J, *invJ, detJ;
558ab72731fSMatthew G. Knepley   PetscReal          xi[4];
559552f7358SJed Brown 
560552f7358SJed Brown   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5649566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
56563a3b9bcSJacob Faibussowitsch   PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5669566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
567ab72731fSMatthew G. Knepley   for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
568ab72731fSMatthew G. Knepley   for (PetscInt d = 0; d < ctx->dim; ++d) {
569552f7358SJed Brown     xi[d] = 0.0;
570ab72731fSMatthew G. Knepley     for (PetscInt f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]);
571ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[(d + 1) * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d];
572552f7358SJed Brown   }
5739566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
5749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5769566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
5773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
578552f7358SJed Brown }
579552f7358SJed Brown 
580ab72731fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v)
581d71ae5a4SJacob Faibussowitsch {
582ab72731fSMatthew G. Knepley   const PetscInt     c        = ctx->cells[p];
583ab72731fSMatthew G. Knepley   const PetscInt     order[3] = {2, 1, 3};
584ab72731fSMatthew G. Knepley   PetscScalar       *x        = NULL;
5857a1931ceSMatthew G. Knepley   PetscReal         *v0, *J, *invJ, detJ;
58656044e6dSMatthew G. Knepley   const PetscScalar *coords;
58756044e6dSMatthew G. Knepley   PetscScalar       *a;
588ab72731fSMatthew G. Knepley   PetscReal          xi[4];
5897a1931ceSMatthew G. Knepley 
5907a1931ceSMatthew G. Knepley   PetscFunctionBegin;
5919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5949566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
59563a3b9bcSJacob Faibussowitsch   PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5969566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
597ab72731fSMatthew G. Knepley   for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
598ab72731fSMatthew G. Knepley   for (PetscInt d = 0; d < ctx->dim; ++d) {
5997a1931ceSMatthew G. Knepley     xi[d] = 0.0;
600ab72731fSMatthew G. Knepley     for (PetscInt f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]);
601ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[order[d] * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d];
6027a1931ceSMatthew G. Knepley   }
6039566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
6049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
6059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
6069566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6087a1931ceSMatthew G. Knepley }
6097a1931ceSMatthew G. Knepley 
610d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
611d71ae5a4SJacob Faibussowitsch {
612552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
613552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
614552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
615552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
616552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
617552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
618552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
619552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
620552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
621552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
622552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
623552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
624552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
625552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
626552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
62756044e6dSMatthew G. Knepley   const PetscScalar *ref;
62856044e6dSMatthew G. Knepley   PetscScalar       *real;
629552f7358SJed Brown 
630552f7358SJed Brown   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
6329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
633552f7358SJed Brown   {
634552f7358SJed Brown     const PetscScalar p0 = ref[0];
635552f7358SJed Brown     const PetscScalar p1 = ref[1];
636552f7358SJed Brown 
637552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
638552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
639552f7358SJed Brown   }
6409566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(28));
6419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
6433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
644552f7358SJed Brown }
645552f7358SJed Brown 
646af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
647d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
648d71ae5a4SJacob Faibussowitsch {
649552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
650552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
651552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
652552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
653552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
654552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
655552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
656552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
657552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
658552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
659552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
66056044e6dSMatthew G. Knepley   const PetscScalar *ref;
661552f7358SJed Brown 
662552f7358SJed Brown   PetscFunctionBegin;
6639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
664552f7358SJed Brown   {
665552f7358SJed Brown     const PetscScalar x       = ref[0];
666552f7358SJed Brown     const PetscScalar y       = ref[1];
667552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
668da80777bSKarl Rupp     PetscScalar       values[4];
669da80777bSKarl Rupp 
6709371c9d4SSatish Balay     values[0] = (x1 - x0 + f_01 * y) * 0.5;
6719371c9d4SSatish Balay     values[1] = (x3 - x0 + f_01 * x) * 0.5;
6729371c9d4SSatish Balay     values[2] = (y1 - y0 + g_01 * y) * 0.5;
6739371c9d4SSatish Balay     values[3] = (y3 - y0 + g_01 * x) * 0.5;
6749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES));
675552f7358SJed Brown   }
6769566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(30));
6779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681552f7358SJed Brown }
682552f7358SJed Brown 
683ab72731fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v)
684d71ae5a4SJacob Faibussowitsch {
685ab72731fSMatthew G. Knepley   const PetscInt     c   = ctx->cells[p];
686de73a395SMatthew G. Knepley   PetscFE            fem = NULL;
687ab72731fSMatthew G. Knepley   DM                 dmCoord;
688552f7358SJed Brown   SNES               snes;
689552f7358SJed Brown   KSP                ksp;
690552f7358SJed Brown   PC                 pc;
691552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
692552f7358SJed Brown   Mat                J;
693716009bfSMatthew G. Knepley   PetscTabulation    T = NULL;
69456044e6dSMatthew G. Knepley   const PetscScalar *coords;
69556044e6dSMatthew G. Knepley   PetscScalar       *a;
696716009bfSMatthew G. Knepley   PetscReal          xir[2] = {0., 0.};
697ab72731fSMatthew G. Knepley   PetscInt           Nf;
6985509d985SMatthew G. Knepley   const PetscInt     dof = ctx->dof;
699ab72731fSMatthew G. Knepley   PetscScalar       *x = NULL, *vertices = NULL;
700ab72731fSMatthew G. Knepley   PetscScalar       *xi;
701ab72731fSMatthew G. Knepley   PetscInt           coordSize, xSize;
702552f7358SJed Brown 
703552f7358SJed Brown   PetscFunctionBegin;
7049566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
705716009bfSMatthew G. Knepley   if (Nf) {
706cfe5744fSMatthew G. Knepley     PetscObject  obj;
707cfe5744fSMatthew G. Knepley     PetscClassId id;
708cfe5744fSMatthew G. Knepley 
7099566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, 0, NULL, &obj));
7109566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
7119371c9d4SSatish Balay     if (id == PETSCFE_CLASSID) {
7129371c9d4SSatish Balay       fem = (PetscFE)obj;
7139371c9d4SSatish Balay       PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T));
7149371c9d4SSatish Balay     }
715716009bfSMatthew G. Knepley   }
7169566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
7179566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
7189566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
7199566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_"));
7209566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
7219566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 2, 2));
7229566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
7239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
7249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
7259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
7269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 2, 2, 2, 2));
7279566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
7289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
7299566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL));
7309566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL));
7319566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7329566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
7339566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
7349566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
735552f7358SJed Brown 
7369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
7379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
7389566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
73963a3b9bcSJacob Faibussowitsch   PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2);
7409566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
7419566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
7429566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
7439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(real, &xi));
744552f7358SJed Brown   xi[0] = coords[p * ctx->dim + 0];
745552f7358SJed Brown   xi[1] = coords[p * ctx->dim + 1];
7469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(real, &xi));
7479566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, real, ref));
7489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ref, &xi));
749cb313848SJed Brown   xir[0] = PetscRealPart(xi[0]);
750cb313848SJed Brown   xir[1] = PetscRealPart(xi[1]);
751cfe5744fSMatthew G. Knepley   if (4 * dof == xSize) {
752ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) * (1 - xir[1]) + x[1 * dof + comp] * xir[0] * (1 - xir[1]) + x[2 * dof + comp] * xir[0] * xir[1] + x[3 * dof + comp] * (1 - xir[0]) * xir[1];
753cfe5744fSMatthew G. Knepley   } else if (dof == xSize) {
754ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
755cfe5744fSMatthew G. Knepley   } else {
7562c71b3e2SJacob Faibussowitsch     PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE");
7579371c9d4SSatish Balay     xir[0] = 2.0 * xir[0] - 1.0;
7589371c9d4SSatish Balay     xir[1] = 2.0 * xir[1] - 1.0;
7599566063dSJacob Faibussowitsch     PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T));
760ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < dof; ++comp) {
7615509d985SMatthew G. Knepley       a[p * dof + comp] = 0.0;
762ab72731fSMatthew G. Knepley       for (PetscInt d = 0; d < xSize / dof; ++d) a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp];
7635509d985SMatthew G. Knepley     }
7645509d985SMatthew G. Knepley   }
7659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ref, &xi));
7669566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
7679566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
7689566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&T));
7699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
7709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
771552f7358SJed Brown 
7729566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
7759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
7769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
7773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
778552f7358SJed Brown }
779552f7358SJed Brown 
780d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
781d71ae5a4SJacob Faibussowitsch {
782552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
783552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
784552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
785552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
7867a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
7877a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
7887a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
789552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
790552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
791552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
7927a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
7937a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
7947a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
795552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
796552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
797552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
798552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
799552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
800552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
801552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
802552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
803552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
804552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
805552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
806552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
807552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
808552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
809552f7358SJed Brown   const PetscScalar  h_1      = z1 - z0;
810552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
811552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
812552f7358SJed Brown   const PetscScalar  h_3      = z3 - z0;
813552f7358SJed Brown   const PetscScalar  f_4      = x4 - x0;
814552f7358SJed Brown   const PetscScalar  g_4      = y4 - y0;
815552f7358SJed Brown   const PetscScalar  h_4      = z4 - z0;
816552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
817552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
818552f7358SJed Brown   const PetscScalar  h_01     = z2 - z1 - z3 + z0;
819552f7358SJed Brown   const PetscScalar  f_12     = x7 - x3 - x4 + x0;
820552f7358SJed Brown   const PetscScalar  g_12     = y7 - y3 - y4 + y0;
821552f7358SJed Brown   const PetscScalar  h_12     = z7 - z3 - z4 + z0;
822552f7358SJed Brown   const PetscScalar  f_02     = x5 - x1 - x4 + x0;
823552f7358SJed Brown   const PetscScalar  g_02     = y5 - y1 - y4 + y0;
824552f7358SJed Brown   const PetscScalar  h_02     = z5 - z1 - z4 + z0;
825552f7358SJed Brown   const PetscScalar  f_012    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
826552f7358SJed Brown   const PetscScalar  g_012    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
827552f7358SJed Brown   const PetscScalar  h_012    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
82856044e6dSMatthew G. Knepley   const PetscScalar *ref;
82956044e6dSMatthew G. Knepley   PetscScalar       *real;
830552f7358SJed Brown 
831552f7358SJed Brown   PetscFunctionBegin;
8329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
8339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
834552f7358SJed Brown   {
835552f7358SJed Brown     const PetscScalar p0 = ref[0];
836552f7358SJed Brown     const PetscScalar p1 = ref[1];
837552f7358SJed Brown     const PetscScalar p2 = ref[2];
838552f7358SJed Brown 
839552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_4 * p2 + f_01 * p0 * p1 + f_12 * p1 * p2 + f_02 * p0 * p2 + f_012 * p0 * p1 * p2;
840552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_4 * p2 + g_01 * p0 * p1 + g_01 * p0 * p1 + g_12 * p1 * p2 + g_02 * p0 * p2 + g_012 * p0 * p1 * p2;
841552f7358SJed Brown     real[2] = z0 + h_1 * p0 + h_3 * p1 + h_4 * p2 + h_01 * p0 * p1 + h_01 * p0 * p1 + h_12 * p1 * p2 + h_02 * p0 * p2 + h_012 * p0 * p1 * p2;
842552f7358SJed Brown   }
8439566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(114));
8449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
8459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
8463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
847552f7358SJed Brown }
848552f7358SJed Brown 
849d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
850d71ae5a4SJacob Faibussowitsch {
851552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
852552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
853552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
854552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
8557a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
8567a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
8577a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
858552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
859552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
860552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
8617a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
8627a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
8637a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
864552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
865552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
866552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
867552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
868552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
869552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
870552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
871552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
872552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
873552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
874552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
875552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
876552f7358SJed Brown   const PetscScalar  f_xy     = x2 - x1 - x3 + x0;
877552f7358SJed Brown   const PetscScalar  g_xy     = y2 - y1 - y3 + y0;
878552f7358SJed Brown   const PetscScalar  h_xy     = z2 - z1 - z3 + z0;
879552f7358SJed Brown   const PetscScalar  f_yz     = x7 - x3 - x4 + x0;
880552f7358SJed Brown   const PetscScalar  g_yz     = y7 - y3 - y4 + y0;
881552f7358SJed Brown   const PetscScalar  h_yz     = z7 - z3 - z4 + z0;
882552f7358SJed Brown   const PetscScalar  f_xz     = x5 - x1 - x4 + x0;
883552f7358SJed Brown   const PetscScalar  g_xz     = y5 - y1 - y4 + y0;
884552f7358SJed Brown   const PetscScalar  h_xz     = z5 - z1 - z4 + z0;
885552f7358SJed Brown   const PetscScalar  f_xyz    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
886552f7358SJed Brown   const PetscScalar  g_xyz    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
887552f7358SJed Brown   const PetscScalar  h_xyz    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
88856044e6dSMatthew G. Knepley   const PetscScalar *ref;
889552f7358SJed Brown 
890552f7358SJed Brown   PetscFunctionBegin;
8919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
892552f7358SJed Brown   {
893552f7358SJed Brown     const PetscScalar x       = ref[0];
894552f7358SJed Brown     const PetscScalar y       = ref[1];
895552f7358SJed Brown     const PetscScalar z       = ref[2];
896552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
897da80777bSKarl Rupp     PetscScalar       values[9];
898da80777bSKarl Rupp 
899da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0;
900da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0;
901da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0;
902da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0;
903da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0;
904da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0;
905da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0;
906da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0;
907da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0;
9081aa26658SKarl Rupp 
9099566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES));
910552f7358SJed Brown   }
9119566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(152));
9129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
9139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
9149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
916552f7358SJed Brown }
917552f7358SJed Brown 
918ab72731fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v)
919d71ae5a4SJacob Faibussowitsch {
920ab72731fSMatthew G. Knepley   const PetscInt     c = ctx->cells[p];
921fafc0619SMatthew G Knepley   DM                 dmCoord;
922552f7358SJed Brown   SNES               snes;
923552f7358SJed Brown   KSP                ksp;
924552f7358SJed Brown   PC                 pc;
925552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
926552f7358SJed Brown   Mat                J;
92756044e6dSMatthew G. Knepley   const PetscScalar *coords;
92856044e6dSMatthew G. Knepley   PetscScalar       *a;
929ab72731fSMatthew G. Knepley   PetscScalar       *x = NULL, *vertices = NULL;
930ab72731fSMatthew G. Knepley   PetscScalar       *xi;
931ab72731fSMatthew G. Knepley   PetscReal          xir[3];
932ab72731fSMatthew G. Knepley   PetscInt           coordSize, xSize;
933552f7358SJed Brown 
934552f7358SJed Brown   PetscFunctionBegin;
9359566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
9369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
9379566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
9389566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_"));
9399566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
9409566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 3, 3));
9419566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
9429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
9439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
9449566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
9459566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 3, 3, 3, 3));
9469566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
9479566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
9489566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL));
9499566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL));
9509566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
9519566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
9529566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
9539566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
954552f7358SJed Brown 
9559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
9569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
9579566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
958cfe5744fSMatthew G. Knepley   PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3);
9599566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
960cfe5744fSMatthew G. Knepley   PetscCheck((8 * ctx->dof == xSize) || (ctx->dof == xSize), ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input closure size %" PetscInt_FMT " should be %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 8 * ctx->dof, ctx->dof);
9619566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
9629566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
9639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(real, &xi));
964552f7358SJed Brown   xi[0] = coords[p * ctx->dim + 0];
965552f7358SJed Brown   xi[1] = coords[p * ctx->dim + 1];
966552f7358SJed Brown   xi[2] = coords[p * ctx->dim + 2];
9679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(real, &xi));
9689566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, real, ref));
9699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ref, &xi));
970cb313848SJed Brown   xir[0] = PetscRealPart(xi[0]);
971cb313848SJed Brown   xir[1] = PetscRealPart(xi[1]);
972cb313848SJed Brown   xir[2] = PetscRealPart(xi[2]);
973cfe5744fSMatthew G. Knepley   if (8 * ctx->dof == xSize) {
974ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < ctx->dof; ++comp) {
9759371c9d4SSatish Balay       a[p * ctx->dof + comp] = x[0 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * (1 - xir[2]) + x[3 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * (1 - xir[2]) + x[2 * ctx->dof + comp] * xir[0] * xir[1] * (1 - xir[2]) + x[1 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * (1 - xir[2]) +
9769371c9d4SSatish Balay                                x[4 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * xir[2] + x[5 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * xir[2] + x[6 * ctx->dof + comp] * xir[0] * xir[1] * xir[2] + x[7 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * xir[2];
977552f7358SJed Brown     }
978cfe5744fSMatthew G. Knepley   } else {
979ab72731fSMatthew G. Knepley     for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
980cfe5744fSMatthew G. Knepley   }
9819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ref, &xi));
9829566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
9839566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
9849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
9859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
986552f7358SJed Brown 
9879566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
9889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
9899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
9909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
9919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
9923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
993552f7358SJed Brown }
994552f7358SJed Brown 
9954267b1a3SMatthew G. Knepley /*@C
9964267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9974267b1a3SMatthew G. Knepley 
998552f7358SJed Brown   Input Parameters:
999f6dfbefdSBarry Smith + ctx - The `DMInterpolationInfo` context
1000f6dfbefdSBarry Smith . dm  - The `DM`
1001552f7358SJed Brown - x   - The local vector containing the field to be interpolated
1002552f7358SJed Brown 
10032fe279fdSBarry Smith   Output Parameter:
1004552f7358SJed Brown . v - The vector containing the interpolated values
10054267b1a3SMatthew G. Knepley 
10064267b1a3SMatthew G. Knepley   Level: beginner
10074267b1a3SMatthew G. Knepley 
10082fe279fdSBarry Smith   Note:
10092fe279fdSBarry Smith   A suitable `v` can be obtained using `DMInterpolationGetVector()`.
10102fe279fdSBarry Smith 
1011f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
10124267b1a3SMatthew G. Knepley @*/
1013d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
1014d71ae5a4SJacob Faibussowitsch {
101566a0a883SMatthew G. Knepley   PetscDS   ds;
101666a0a883SMatthew G. Knepley   PetscInt  n, p, Nf, field;
101766a0a883SMatthew G. Knepley   PetscBool useDS = PETSC_FALSE;
1018552f7358SJed Brown 
1019552f7358SJed Brown   PetscFunctionBegin;
1020552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1021552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1022552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
10239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
102463a3b9bcSJacob Faibussowitsch   PetscCheck(n == ctx->n * ctx->dof, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %" PetscInt_FMT " should be %" PetscInt_FMT, n, ctx->n * ctx->dof);
10253ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
10269566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1027680d18d5SMatthew G. Knepley   if (ds) {
102866a0a883SMatthew G. Knepley     useDS = PETSC_TRUE;
10299566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &Nf));
1030680d18d5SMatthew G. Knepley     for (field = 0; field < Nf; ++field) {
103166a0a883SMatthew G. Knepley       PetscObject  obj;
1032680d18d5SMatthew G. Knepley       PetscClassId id;
1033680d18d5SMatthew G. Knepley 
10349566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, field, &obj));
10359566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
1036*c3af198dSMatthew G. Knepley       if (id != PETSCFE_CLASSID && id != PETSCFV_CLASSID) {
10379371c9d4SSatish Balay         useDS = PETSC_FALSE;
10389371c9d4SSatish Balay         break;
10399371c9d4SSatish Balay       }
104066a0a883SMatthew G. Knepley     }
104166a0a883SMatthew G. Knepley   }
104266a0a883SMatthew G. Knepley   if (useDS) {
104366a0a883SMatthew G. Knepley     const PetscScalar *coords;
104466a0a883SMatthew G. Knepley     PetscScalar       *interpolant;
104566a0a883SMatthew G. Knepley     PetscInt           cdim, d;
104666a0a883SMatthew G. Knepley 
10479566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm, &cdim));
10489566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ctx->coords, &coords));
10499566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(v, &interpolant));
1050680d18d5SMatthew G. Knepley     for (p = 0; p < ctx->n; ++p) {
105166a0a883SMatthew G. Knepley       PetscReal    pcoords[3], xi[3];
1052680d18d5SMatthew G. Knepley       PetscScalar *xa   = NULL;
105366a0a883SMatthew G. Knepley       PetscInt     coff = 0, foff = 0, clSize;
1054680d18d5SMatthew G. Knepley 
105552aa1562SMatthew G. Knepley       if (ctx->cells[p] < 0) continue;
1056680d18d5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]);
10579566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi));
10589566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
105966a0a883SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
106066a0a883SMatthew G. Knepley         PetscTabulation T;
1061*c3af198dSMatthew G. Knepley         PetscObject     obj;
1062*c3af198dSMatthew G. Knepley         PetscClassId    id;
106366a0a883SMatthew G. Knepley 
1064*c3af198dSMatthew G. Knepley         PetscCall(PetscDSGetDiscretization(ds, field, &obj));
1065*c3af198dSMatthew G. Knepley         PetscCall(PetscObjectGetClassId(obj, &id));
1066*c3af198dSMatthew G. Knepley         if (id == PETSCFE_CLASSID) {
1067*c3af198dSMatthew G. Knepley           PetscFE fe = (PetscFE)obj;
1068*c3af198dSMatthew G. Knepley 
10699566063dSJacob Faibussowitsch           PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T));
1070680d18d5SMatthew G. Knepley           {
1071680d18d5SMatthew G. Knepley             const PetscReal *basis = T->T[0];
1072680d18d5SMatthew G. Knepley             const PetscInt   Nb    = T->Nb;
1073680d18d5SMatthew G. Knepley             const PetscInt   Nc    = T->Nc;
107466a0a883SMatthew G. Knepley 
1075*c3af198dSMatthew G. Knepley             for (PetscInt fc = 0; fc < Nc; ++fc) {
107666a0a883SMatthew G. Knepley               interpolant[p * ctx->dof + coff + fc] = 0.0;
1077*c3af198dSMatthew G. Knepley               for (PetscInt f = 0; f < Nb; ++f) interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc];
1078680d18d5SMatthew G. Knepley             }
107966a0a883SMatthew G. Knepley             coff += Nc;
108066a0a883SMatthew G. Knepley             foff += Nb;
1081680d18d5SMatthew G. Knepley           }
10829566063dSJacob Faibussowitsch           PetscCall(PetscTabulationDestroy(&T));
1083*c3af198dSMatthew G. Knepley         } else if (id == PETSCFV_CLASSID) {
1084*c3af198dSMatthew G. Knepley           PetscFV  fv = (PetscFV)obj;
1085*c3af198dSMatthew G. Knepley           PetscInt Nc;
1086*c3af198dSMatthew G. Knepley 
1087*c3af198dSMatthew G. Knepley           // TODO Could use reconstruction if available
1088*c3af198dSMatthew G. Knepley           PetscCall(PetscFVGetNumComponents(fv, &Nc));
1089*c3af198dSMatthew G. Knepley           for (PetscInt fc = 0; fc < Nc; ++fc) interpolant[p * ctx->dof + coff + fc] = xa[foff + fc];
1090*c3af198dSMatthew G. Knepley           coff += Nc;
1091*c3af198dSMatthew G. Knepley           foff += Nc;
1092*c3af198dSMatthew G. Knepley         }
1093680d18d5SMatthew G. Knepley       }
10949566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
109563a3b9bcSJacob Faibussowitsch       PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof);
1096*c3af198dSMatthew G. Knepley       PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE/FV space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize);
109766a0a883SMatthew G. Knepley     }
10989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
10999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(v, &interpolant));
110066a0a883SMatthew G. Knepley   } else {
1101ab72731fSMatthew G. Knepley     for (PetscInt p = 0; p < ctx->n; ++p) {
1102ab72731fSMatthew G. Knepley       const PetscInt cell = ctx->cells[p];
110366a0a883SMatthew G. Knepley       DMPolytopeType ct;
110466a0a883SMatthew G. Knepley 
1105ab72731fSMatthew G. Knepley       PetscCall(DMPlexGetCellType(dm, cell, &ct));
1106680d18d5SMatthew G. Knepley       switch (ct) {
1107d71ae5a4SJacob Faibussowitsch       case DM_POLYTOPE_SEGMENT:
1108ab72731fSMatthew G. Knepley         PetscCall(DMInterpolate_Segment_Private(ctx, dm, p, x, v));
1109d71ae5a4SJacob Faibussowitsch         break;
1110d71ae5a4SJacob Faibussowitsch       case DM_POLYTOPE_TRIANGLE:
1111ab72731fSMatthew G. Knepley         PetscCall(DMInterpolate_Triangle_Private(ctx, dm, p, x, v));
1112d71ae5a4SJacob Faibussowitsch         break;
1113d71ae5a4SJacob Faibussowitsch       case DM_POLYTOPE_QUADRILATERAL:
1114ab72731fSMatthew G. Knepley         PetscCall(DMInterpolate_Quad_Private(ctx, dm, p, x, v));
1115d71ae5a4SJacob Faibussowitsch         break;
1116d71ae5a4SJacob Faibussowitsch       case DM_POLYTOPE_TETRAHEDRON:
1117ab72731fSMatthew G. Knepley         PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, p, x, v));
1118d71ae5a4SJacob Faibussowitsch         break;
1119d71ae5a4SJacob Faibussowitsch       case DM_POLYTOPE_HEXAHEDRON:
1120ab72731fSMatthew G. Knepley         PetscCall(DMInterpolate_Hex_Private(ctx, dm, cell, x, v));
1121d71ae5a4SJacob Faibussowitsch         break;
1122d71ae5a4SJacob Faibussowitsch       default:
1123d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]);
1124680d18d5SMatthew G. Knepley       }
1125552f7358SJed Brown     }
1126ab72731fSMatthew G. Knepley   }
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1128552f7358SJed Brown }
1129552f7358SJed Brown 
11304267b1a3SMatthew G. Knepley /*@C
1131f6dfbefdSBarry Smith   DMInterpolationDestroy - Destroys a `DMInterpolationInfo` context
11324267b1a3SMatthew G. Knepley 
1133c3339decSBarry Smith   Collective
11344267b1a3SMatthew G. Knepley 
11354267b1a3SMatthew G. Knepley   Input Parameter:
11364267b1a3SMatthew G. Knepley . ctx - the context
11374267b1a3SMatthew G. Knepley 
11384267b1a3SMatthew G. Knepley   Level: beginner
11394267b1a3SMatthew G. Knepley 
1140f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
11414267b1a3SMatthew G. Knepley @*/
1142d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
1143d71ae5a4SJacob Faibussowitsch {
1144552f7358SJed Brown   PetscFunctionBegin;
11454f572ea9SToby Isaac   PetscAssertPointer(ctx, 1);
11469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->coords));
11479566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->points));
11489566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->cells));
11499566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
11500298fd71SBarry Smith   *ctx = NULL;
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152552f7358SJed Brown }
1153cc0c4584SMatthew G. Knepley 
1154cc0c4584SMatthew G. Knepley /*@C
1155cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1156cc0c4584SMatthew G. Knepley 
1157c3339decSBarry Smith   Collective
1158cc0c4584SMatthew G. Knepley 
1159cc0c4584SMatthew G. Knepley   Input Parameters:
1160f6dfbefdSBarry Smith + snes   - the `SNES` context
1161cc0c4584SMatthew G. Knepley . its    - iteration number
1162cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1163f6dfbefdSBarry Smith - vf     - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII`
1164cc0c4584SMatthew G. Knepley 
11652fe279fdSBarry Smith   Level: intermediate
11662fe279fdSBarry Smith 
1167f6dfbefdSBarry Smith   Note:
1168cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1169cc0c4584SMatthew G. Knepley 
1170f6dfbefdSBarry Smith .seealso: `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()`
1171cc0c4584SMatthew G. Knepley @*/
1172d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1173d71ae5a4SJacob Faibussowitsch {
1174d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1175cc0c4584SMatthew G. Knepley   Vec                res;
1176cc0c4584SMatthew G. Knepley   DM                 dm;
1177cc0c4584SMatthew G. Knepley   PetscSection       s;
1178cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1179cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1180cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1181cc0c4584SMatthew G. Knepley 
1182cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11834d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
11849566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
11859566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
11869566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
11879566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
11889566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
11899566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms));
11909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(res, &r));
1191cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1192cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1193cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1194cc0c4584SMatthew G. Knepley 
11959566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
11969566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
1197cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d]));
1198cc0c4584SMatthew G. Knepley     }
1199cc0c4584SMatthew G. Knepley   }
12009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(res, &r));
12011c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
12029566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
12039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
120463a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm));
1205cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
12069566063dSJacob Faibussowitsch     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
12079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f])));
1208cc0c4584SMatthew G. Knepley   }
12099566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
12109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
12119566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
12129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lnorms, norms));
12133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1214cc0c4584SMatthew G. Knepley }
121524cdb843SMatthew G. Knepley 
121624cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
121724cdb843SMatthew G. Knepley 
121824cdb843SMatthew G. Knepley /*@
12198db23a53SJed Brown   DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user
122024cdb843SMatthew G. Knepley 
122124cdb843SMatthew G. Knepley   Input Parameters:
122224cdb843SMatthew G. Knepley + dm   - The mesh
122324cdb843SMatthew G. Knepley . X    - Local solution
122424cdb843SMatthew G. Knepley - user - The user context
122524cdb843SMatthew G. Knepley 
122624cdb843SMatthew G. Knepley   Output Parameter:
122724cdb843SMatthew G. Knepley . F - Local output vector
122824cdb843SMatthew G. Knepley 
12292fe279fdSBarry Smith   Level: developer
12302fe279fdSBarry Smith 
1231f6dfbefdSBarry Smith   Note:
1232f6dfbefdSBarry Smith   The residual is summed into F; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in F is intentional.
12338db23a53SJed Brown 
1234f6dfbefdSBarry Smith .seealso: `DM`, `DMPlexComputeJacobianAction()`
123524cdb843SMatthew G. Knepley @*/
1236d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1237d71ae5a4SJacob Faibussowitsch {
12386da023fcSToby Isaac   DM       plex;
1239083401c6SMatthew G. Knepley   IS       allcellIS;
12406528b96dSMatthew G. Knepley   PetscInt Nds, s;
124124cdb843SMatthew G. Knepley 
124224cdb843SMatthew G. Knepley   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12459566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
12466528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
12476528b96dSMatthew G. Knepley     PetscDS      ds;
12486528b96dSMatthew G. Knepley     IS           cellIS;
124906ad1575SMatthew G. Knepley     PetscFormKey key;
12506528b96dSMatthew G. Knepley 
125107218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
12526528b96dSMatthew G. Knepley     key.value = 0;
12536528b96dSMatthew G. Knepley     key.field = 0;
125406ad1575SMatthew G. Knepley     key.part  = 0;
12556528b96dSMatthew G. Knepley     if (!key.label) {
12569566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
12576528b96dSMatthew G. Knepley       cellIS = allcellIS;
12586528b96dSMatthew G. Knepley     } else {
12596528b96dSMatthew G. Knepley       IS pointIS;
12606528b96dSMatthew G. Knepley 
12616528b96dSMatthew G. Knepley       key.value = 1;
12629566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
12639566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12649566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
12656528b96dSMatthew G. Knepley     }
12669566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
12679566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
12686528b96dSMatthew G. Knepley   }
12699566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
12709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
12713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12726528b96dSMatthew G. Knepley }
12736528b96dSMatthew G. Knepley 
1274bdd6f66aSToby Isaac /*@
1275d2b2dc1eSMatthew G. Knepley   DMPlexSNESComputeResidualDS - Sums the local residual into vector F from the local input X using all pointwise functions with unique keys in the DS
1276d2b2dc1eSMatthew G. Knepley 
1277d2b2dc1eSMatthew G. Knepley   Input Parameters:
1278d2b2dc1eSMatthew G. Knepley + dm   - The mesh
1279d2b2dc1eSMatthew G. Knepley . X    - Local solution
1280d2b2dc1eSMatthew G. Knepley - user - The user context
1281d2b2dc1eSMatthew G. Knepley 
1282d2b2dc1eSMatthew G. Knepley   Output Parameter:
1283d2b2dc1eSMatthew G. Knepley . F - Local output vector
1284d2b2dc1eSMatthew G. Knepley 
1285d2b2dc1eSMatthew G. Knepley   Level: developer
1286d2b2dc1eSMatthew G. Knepley 
1287d2b2dc1eSMatthew G. Knepley   Note:
1288d2b2dc1eSMatthew G. Knepley   The residual is summed into F; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in F is intentional.
1289d2b2dc1eSMatthew G. Knepley 
1290d2b2dc1eSMatthew G. Knepley .seealso: `DM`, `DMPlexComputeJacobianAction()`
1291d2b2dc1eSMatthew G. Knepley @*/
1292d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualDS(DM dm, Vec X, Vec F, void *user)
1293d2b2dc1eSMatthew G. Knepley {
1294d2b2dc1eSMatthew G. Knepley   DM       plex;
1295d2b2dc1eSMatthew G. Knepley   IS       allcellIS;
1296d2b2dc1eSMatthew G. Knepley   PetscInt Nds, s;
1297d2b2dc1eSMatthew G. Knepley 
1298d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
1299d2b2dc1eSMatthew G. Knepley   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
1300d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1301d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetNumDS(dm, &Nds));
1302d2b2dc1eSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1303d2b2dc1eSMatthew G. Knepley     PetscDS ds;
1304d2b2dc1eSMatthew G. Knepley     DMLabel label;
1305d2b2dc1eSMatthew G. Knepley     IS      cellIS;
1306d2b2dc1eSMatthew G. Knepley 
1307d2b2dc1eSMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL));
1308d2b2dc1eSMatthew G. Knepley     {
1309d2b2dc1eSMatthew G. Knepley       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
1310d2b2dc1eSMatthew G. Knepley       PetscWeakForm     wf;
1311d2b2dc1eSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
1312d2b2dc1eSMatthew G. Knepley       PetscFormKey     *reskeys;
1313d2b2dc1eSMatthew G. Knepley 
1314d2b2dc1eSMatthew G. Knepley       /* Get unique residual keys */
1315d2b2dc1eSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
1316d2b2dc1eSMatthew G. Knepley         PetscInt Nkm;
1317d2b2dc1eSMatthew G. Knepley         PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm));
1318d2b2dc1eSMatthew G. Knepley         Nk += Nkm;
1319d2b2dc1eSMatthew G. Knepley       }
1320d2b2dc1eSMatthew G. Knepley       PetscCall(PetscMalloc1(Nk, &reskeys));
1321d2b2dc1eSMatthew G. Knepley       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys));
1322d2b2dc1eSMatthew G. Knepley       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
1323d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFormKeySort(Nk, reskeys));
1324d2b2dc1eSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
1325d2b2dc1eSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
1326d2b2dc1eSMatthew G. Knepley           ++k;
1327d2b2dc1eSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
1328d2b2dc1eSMatthew G. Knepley         }
1329d2b2dc1eSMatthew G. Knepley       }
1330d2b2dc1eSMatthew G. Knepley       Nk = k;
1331d2b2dc1eSMatthew G. Knepley 
1332d2b2dc1eSMatthew G. Knepley       PetscCall(PetscDSGetWeakForm(ds, &wf));
1333d2b2dc1eSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
1334d2b2dc1eSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
1335d2b2dc1eSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
1336d2b2dc1eSMatthew G. Knepley 
1337d2b2dc1eSMatthew G. Knepley         if (!label) {
1338d2b2dc1eSMatthew G. Knepley           PetscCall(PetscObjectReference((PetscObject)allcellIS));
1339d2b2dc1eSMatthew G. Knepley           cellIS = allcellIS;
1340d2b2dc1eSMatthew G. Knepley         } else {
1341d2b2dc1eSMatthew G. Knepley           IS pointIS;
1342d2b2dc1eSMatthew G. Knepley 
1343d2b2dc1eSMatthew G. Knepley           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
1344d2b2dc1eSMatthew G. Knepley           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
1345d2b2dc1eSMatthew G. Knepley           PetscCall(ISDestroy(&pointIS));
1346d2b2dc1eSMatthew G. Knepley         }
1347d2b2dc1eSMatthew G. Knepley         PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
1348d2b2dc1eSMatthew G. Knepley         PetscCall(ISDestroy(&cellIS));
1349d2b2dc1eSMatthew G. Knepley       }
1350d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFree(reskeys));
1351d2b2dc1eSMatthew G. Knepley     }
1352d2b2dc1eSMatthew G. Knepley   }
1353d2b2dc1eSMatthew G. Knepley   PetscCall(ISDestroy(&allcellIS));
1354d2b2dc1eSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
1355d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1356d2b2dc1eSMatthew G. Knepley }
1357d2b2dc1eSMatthew G. Knepley 
1358d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
1359d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualCEED(DM dm, Vec locX, Vec locF, void *user)
1360d2b2dc1eSMatthew G. Knepley {
1361d2b2dc1eSMatthew G. Knepley   Ceed       ceed;
1362d2b2dc1eSMatthew G. Knepley   DMCeed     sd = dm->dmceed;
1363d2b2dc1eSMatthew G. Knepley   CeedVector clocX, clocF;
1364d2b2dc1eSMatthew G. Knepley 
1365d2b2dc1eSMatthew G. Knepley   PetscFunctionBegin;
1366d2b2dc1eSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
1367d2b2dc1eSMatthew G. Knepley   PetscCheck(sd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DM has no CEED data. Call DMCeedCreate() before computing the residual.");
1368d2b2dc1eSMatthew G. Knepley   PetscCall(DMCeedComputeGeometry(dm, sd));
1369d2b2dc1eSMatthew G. Knepley 
1370d2b2dc1eSMatthew G. Knepley   PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX));
1371d2b2dc1eSMatthew G. Knepley   PetscCall(VecGetCeedVector(locF, ceed, &clocF));
1372d2b2dc1eSMatthew G. Knepley   PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE));
1373d2b2dc1eSMatthew G. Knepley   PetscCall(VecRestoreCeedVectorRead(locX, &clocX));
1374d2b2dc1eSMatthew G. Knepley   PetscCall(VecRestoreCeedVector(locF, &clocF));
1375d2b2dc1eSMatthew G. Knepley 
1376d2b2dc1eSMatthew G. Knepley   {
1377d2b2dc1eSMatthew G. Knepley     DM_Plex *mesh = (DM_Plex *)dm->data;
1378d2b2dc1eSMatthew G. Knepley 
1379d2b2dc1eSMatthew G. Knepley     if (mesh->printFEM) {
1380d2b2dc1eSMatthew G. Knepley       PetscSection section;
1381d2b2dc1eSMatthew G. Knepley       Vec          locFbc;
1382d2b2dc1eSMatthew G. Knepley       PetscInt     pStart, pEnd, p, maxDof;
1383d2b2dc1eSMatthew G. Knepley       PetscScalar *zeroes;
1384d2b2dc1eSMatthew G. Knepley 
1385d2b2dc1eSMatthew G. Knepley       PetscCall(DMGetLocalSection(dm, &section));
1386d2b2dc1eSMatthew G. Knepley       PetscCall(VecDuplicate(locF, &locFbc));
1387d2b2dc1eSMatthew G. Knepley       PetscCall(VecCopy(locF, locFbc));
1388d2b2dc1eSMatthew G. Knepley       PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1389d2b2dc1eSMatthew G. Knepley       PetscCall(PetscSectionGetMaxDof(section, &maxDof));
1390d2b2dc1eSMatthew G. Knepley       PetscCall(PetscCalloc1(maxDof, &zeroes));
1391d2b2dc1eSMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
1392d2b2dc1eSMatthew G. Knepley       PetscCall(PetscFree(zeroes));
1393d2b2dc1eSMatthew G. Knepley       PetscCall(DMPrintLocalVec(dm, "Residual", mesh->printTol, locFbc));
1394d2b2dc1eSMatthew G. Knepley       PetscCall(VecDestroy(&locFbc));
1395d2b2dc1eSMatthew G. Knepley     }
1396d2b2dc1eSMatthew G. Knepley   }
1397d2b2dc1eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1398d2b2dc1eSMatthew G. Knepley }
1399d2b2dc1eSMatthew G. Knepley #endif
1400d2b2dc1eSMatthew G. Knepley 
1401d2b2dc1eSMatthew G. Knepley /*@
1402bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1403bdd6f66aSToby Isaac 
1404bdd6f66aSToby Isaac   Input Parameters:
1405bdd6f66aSToby Isaac + dm   - The mesh
1406bdd6f66aSToby Isaac - user - The user context
1407bdd6f66aSToby Isaac 
1408bdd6f66aSToby Isaac   Output Parameter:
1409bdd6f66aSToby Isaac . X - Local solution
1410bdd6f66aSToby Isaac 
1411bdd6f66aSToby Isaac   Level: developer
1412bdd6f66aSToby Isaac 
1413f6dfbefdSBarry Smith .seealso: `DMPLEX`, `DMPlexComputeJacobianAction()`
1414bdd6f66aSToby Isaac @*/
1415d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1416d71ae5a4SJacob Faibussowitsch {
1417bdd6f66aSToby Isaac   DM plex;
1418bdd6f66aSToby Isaac 
1419bdd6f66aSToby Isaac   PetscFunctionBegin;
14209566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
14219566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL));
14229566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
14233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1424bdd6f66aSToby Isaac }
1425bdd6f66aSToby Isaac 
14267a73cf09SMatthew G. Knepley /*@
14278e3a2eefSMatthew G. Knepley   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
14287a73cf09SMatthew G. Knepley 
14297a73cf09SMatthew G. Knepley   Input Parameters:
1430f6dfbefdSBarry Smith + dm   - The `DM`
14317a73cf09SMatthew G. Knepley . X    - Local solution vector
14327a73cf09SMatthew G. Knepley . Y    - Local input vector
14337a73cf09SMatthew G. Knepley - user - The user context
14347a73cf09SMatthew G. Knepley 
14357a73cf09SMatthew G. Knepley   Output Parameter:
143669d47153SPierre Jolivet . F - local output vector
14377a73cf09SMatthew G. Knepley 
14387a73cf09SMatthew G. Knepley   Level: developer
14397a73cf09SMatthew G. Knepley 
14408e3a2eefSMatthew G. Knepley   Notes:
1441f6dfbefdSBarry Smith   Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly.
14428e3a2eefSMatthew G. Knepley 
1443f6dfbefdSBarry Smith .seealso: `DM`, ``DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()`
14447a73cf09SMatthew G. Knepley @*/
1445d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
1446d71ae5a4SJacob Faibussowitsch {
14478e3a2eefSMatthew G. Knepley   DM       plex;
14488e3a2eefSMatthew G. Knepley   IS       allcellIS;
14498e3a2eefSMatthew G. Knepley   PetscInt Nds, s;
1450a925c78cSMatthew G. Knepley 
1451a925c78cSMatthew G. Knepley   PetscFunctionBegin;
14529566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
14539566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
14549566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
14558e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
14568e3a2eefSMatthew G. Knepley     PetscDS ds;
14578e3a2eefSMatthew G. Knepley     DMLabel label;
14588e3a2eefSMatthew G. Knepley     IS      cellIS;
14597a73cf09SMatthew G. Knepley 
146007218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL));
14618e3a2eefSMatthew G. Knepley     {
146206ad1575SMatthew G. Knepley       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
14638e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
14648e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
146506ad1575SMatthew G. Knepley       PetscFormKey     *jackeys;
14668e3a2eefSMatthew G. Knepley 
14678e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
14688e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
14698e3a2eefSMatthew G. Knepley         PetscInt Nkm;
14709566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm));
14718e3a2eefSMatthew G. Knepley         Nk += Nkm;
14728e3a2eefSMatthew G. Knepley       }
14739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &jackeys));
147448a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys));
147563a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
14769566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, jackeys));
14778e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
14788e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
14798e3a2eefSMatthew G. Knepley           ++k;
14808e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
14818e3a2eefSMatthew G. Knepley         }
14828e3a2eefSMatthew G. Knepley       }
14838e3a2eefSMatthew G. Knepley       Nk = k;
14848e3a2eefSMatthew G. Knepley 
14859566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
14868e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
14878e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
14888e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
14898e3a2eefSMatthew G. Knepley 
14908e3a2eefSMatthew G. Knepley         if (!label) {
14919566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
14928e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
14937a73cf09SMatthew G. Knepley         } else {
14948e3a2eefSMatthew G. Knepley           IS pointIS;
1495a925c78cSMatthew G. Knepley 
14969566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
14979566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14989566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
1499a925c78cSMatthew G. Knepley         }
15009566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user));
15019566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
15028e3a2eefSMatthew G. Knepley       }
15039566063dSJacob Faibussowitsch       PetscCall(PetscFree(jackeys));
15048e3a2eefSMatthew G. Knepley     }
15058e3a2eefSMatthew G. Knepley   }
15069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
15079566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
15083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1509a925c78cSMatthew G. Knepley }
1510a925c78cSMatthew G. Knepley 
151124cdb843SMatthew G. Knepley /*@
15122fe279fdSBarry Smith   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix `Jac` at the local solution `X` using pointwise functions specified by the user.
151324cdb843SMatthew G. Knepley 
151424cdb843SMatthew G. Knepley   Input Parameters:
15152fe279fdSBarry Smith + dm   - The `DM`
151624cdb843SMatthew G. Knepley . X    - Local input vector
151724cdb843SMatthew G. Knepley - user - The user context
151824cdb843SMatthew G. Knepley 
15192fe279fdSBarry Smith   Output Parameters:
15202fe279fdSBarry Smith + Jac  - Jacobian matrix
15212fe279fdSBarry Smith - JacP - approximate Jacobian from which the preconditioner will be built, often `Jac`
15222fe279fdSBarry Smith 
15232fe279fdSBarry Smith   Level: developer
152424cdb843SMatthew G. Knepley 
152524cdb843SMatthew G. Knepley   Note:
152624cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
152724cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
152824cdb843SMatthew G. Knepley 
1529f6dfbefdSBarry Smith .seealso: `DMPLEX`, `Mat`
153024cdb843SMatthew G. Knepley @*/
1531d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user)
1532d71ae5a4SJacob Faibussowitsch {
15336da023fcSToby Isaac   DM        plex;
1534083401c6SMatthew G. Knepley   IS        allcellIS;
1535f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
15366528b96dSMatthew G. Knepley   PetscInt  Nds, s;
153724cdb843SMatthew G. Knepley 
153824cdb843SMatthew G. Knepley   PetscFunctionBegin;
15399566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
15409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
15419566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1542083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1543083401c6SMatthew G. Knepley     PetscDS      ds;
1544083401c6SMatthew G. Knepley     IS           cellIS;
154506ad1575SMatthew G. Knepley     PetscFormKey key;
1546083401c6SMatthew G. Knepley 
154707218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
15486528b96dSMatthew G. Knepley     key.value = 0;
15496528b96dSMatthew G. Knepley     key.field = 0;
155006ad1575SMatthew G. Knepley     key.part  = 0;
15516528b96dSMatthew G. Knepley     if (!key.label) {
15529566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1553083401c6SMatthew G. Knepley       cellIS = allcellIS;
1554083401c6SMatthew G. Knepley     } else {
1555083401c6SMatthew G. Knepley       IS pointIS;
1556083401c6SMatthew G. Knepley 
15576528b96dSMatthew G. Knepley       key.value = 1;
15589566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
15599566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
15609566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1561083401c6SMatthew G. Knepley     }
1562083401c6SMatthew G. Knepley     if (!s) {
15639566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
15649566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
15659566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
15669566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
1567083401c6SMatthew G. Knepley     }
15689566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user));
15699566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1570083401c6SMatthew G. Knepley   }
15719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
15729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
15733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157424cdb843SMatthew G. Knepley }
15751878804aSMatthew G. Knepley 
15769371c9d4SSatish Balay struct _DMSNESJacobianMFCtx {
15778e3a2eefSMatthew G. Knepley   DM    dm;
15788e3a2eefSMatthew G. Knepley   Vec   X;
15798e3a2eefSMatthew G. Knepley   void *ctx;
15808e3a2eefSMatthew G. Knepley };
15818e3a2eefSMatthew G. Knepley 
1582d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
1583d71ae5a4SJacob Faibussowitsch {
15848e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15858e3a2eefSMatthew G. Knepley 
15868e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15879566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
15889566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, NULL));
15899566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
15909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->X));
15919566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
15923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15938e3a2eefSMatthew G. Knepley }
15948e3a2eefSMatthew G. Knepley 
1595d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
1596d71ae5a4SJacob Faibussowitsch {
15978e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15988e3a2eefSMatthew G. Knepley 
15998e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
16009566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
16019566063dSJacob Faibussowitsch   PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx));
16023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16038e3a2eefSMatthew G. Knepley }
16048e3a2eefSMatthew G. Knepley 
16058e3a2eefSMatthew G. Knepley /*@
1606f6dfbefdSBarry Smith   DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free
16078e3a2eefSMatthew G. Knepley 
1608c3339decSBarry Smith   Collective
16098e3a2eefSMatthew G. Knepley 
16108e3a2eefSMatthew G. Knepley   Input Parameters:
1611f6dfbefdSBarry Smith + dm   - The `DM`
16128e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
16132fe279fdSBarry Smith - user - A user context, or `NULL`
16148e3a2eefSMatthew G. Knepley 
16158e3a2eefSMatthew G. Knepley   Output Parameter:
1616f6dfbefdSBarry Smith . J - The `Mat`
16178e3a2eefSMatthew G. Knepley 
16188e3a2eefSMatthew G. Knepley   Level: advanced
16198e3a2eefSMatthew G. Knepley 
1620f6dfbefdSBarry Smith   Note:
16212fe279fdSBarry Smith   Vec `X` is kept in `J`, so updating `X` then updates the evaluation point.
16228e3a2eefSMatthew G. Knepley 
1623f6dfbefdSBarry Smith .seealso: `DM`, `DMSNESComputeJacobianAction()`
16248e3a2eefSMatthew G. Knepley @*/
1625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
1626d71ae5a4SJacob Faibussowitsch {
16278e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
16288e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
16298e3a2eefSMatthew G. Knepley 
16308e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
16319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
16329566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
16339566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
16349566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
16359566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, n, n, N, N));
16369566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
16379566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)X));
16389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &ctx));
16398e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
16408e3a2eefSMatthew G. Knepley   ctx->X   = X;
16418e3a2eefSMatthew G. Knepley   ctx->ctx = user;
16429566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, ctx));
16439566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private));
16449566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private));
16453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16468e3a2eefSMatthew G. Knepley }
16478e3a2eefSMatthew G. Knepley 
164838cfc46eSPierre Jolivet /*
164938cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
165038cfc46eSPierre Jolivet 
165138cfc46eSPierre Jolivet    Input Parameters:
16522fe279fdSBarry Smith +     X - `SNES` linearization point
165338cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
165438cfc46eSPierre Jolivet 
165538cfc46eSPierre Jolivet    Output Parameter:
165638cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
165738cfc46eSPierre Jolivet 
165838cfc46eSPierre Jolivet    Level: intermediate
165938cfc46eSPierre Jolivet 
1660db781477SPatrick Sanan .seealso: `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()`
166138cfc46eSPierre Jolivet */
1662d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
1663d71ae5a4SJacob Faibussowitsch {
166438cfc46eSPierre Jolivet   SNES   snes;
166538cfc46eSPierre Jolivet   Mat    pJ;
166638cfc46eSPierre Jolivet   DM     ovldm, origdm;
166738cfc46eSPierre Jolivet   DMSNES sdm;
166838cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM, Vec, void *);
166938cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *);
167038cfc46eSPierre Jolivet   void *bctx, *jctx;
167138cfc46eSPierre Jolivet 
167238cfc46eSPierre Jolivet   PetscFunctionBegin;
16739566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ));
167428b400f6SJacob Faibussowitsch   PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat");
16759566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm));
167628b400f6SJacob Faibussowitsch   PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM");
16779566063dSJacob Faibussowitsch   PetscCall(MatGetDM(pJ, &ovldm));
16789566063dSJacob Faibussowitsch   PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx));
16799566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx));
16809566063dSJacob Faibussowitsch   PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx));
16819566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx));
16829566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes));
168338cfc46eSPierre Jolivet   if (!snes) {
16849566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes));
16859566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, ovldm));
16869566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes));
16879566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)snes));
168838cfc46eSPierre Jolivet   }
16899566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(ovldm, &sdm));
16909566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(X));
1691800f99ffSJeremy L Thompson   {
1692800f99ffSJeremy L Thompson     void *ctx;
1693800f99ffSJeremy L Thompson     PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *);
1694800f99ffSJeremy L Thompson     PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx));
1695800f99ffSJeremy L Thompson     PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx));
1696800f99ffSJeremy L Thompson   }
16979566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(X));
169838cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
169938cfc46eSPierre Jolivet   {
170038cfc46eSPierre Jolivet     Mat locpJ;
170138cfc46eSPierre Jolivet 
17029566063dSJacob Faibussowitsch     PetscCall(MatISGetLocalMat(pJ, &locpJ));
17039566063dSJacob Faibussowitsch     PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN));
170438cfc46eSPierre Jolivet   }
17053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170638cfc46eSPierre Jolivet }
170738cfc46eSPierre Jolivet 
1708a925c78cSMatthew G. Knepley /*@
1709f6dfbefdSBarry Smith   DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, residual, and Jacobian.
17109f520fc2SToby Isaac 
17119f520fc2SToby Isaac   Input Parameters:
1712f6dfbefdSBarry Smith + dm          - The `DM` object
1713f6dfbefdSBarry Smith . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see `PetscDSAddBoundary()`)
1714f6dfbefdSBarry Smith . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see `PetscDSSetResidual()`)
1715f6dfbefdSBarry Smith - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see `PetscDSSetJacobian()`)
17161a244344SSatish Balay 
17171a244344SSatish Balay   Level: developer
1718f6dfbefdSBarry Smith 
1719f6dfbefdSBarry Smith .seealso: `DMPLEX`, `SNES`
17209f520fc2SToby Isaac @*/
1721d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
1722d71ae5a4SJacob Faibussowitsch {
1723d2b2dc1eSMatthew G. Knepley   PetscBool useCeed;
1724d2b2dc1eSMatthew G. Knepley 
17259f520fc2SToby Isaac   PetscFunctionBegin;
1726d2b2dc1eSMatthew G. Knepley   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
17279566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, boundaryctx));
1728d2b2dc1eSMatthew G. Knepley   if (useCeed) {
1729d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
1730d2b2dc1eSMatthew G. Knepley     PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualCEED, residualctx));
1731d2b2dc1eSMatthew G. Knepley #else
1732d2b2dc1eSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot use CEED traversals without LibCEED. Rerun configure with --download-ceed");
1733d2b2dc1eSMatthew G. Knepley #endif
1734d2b2dc1eSMatthew G. Knepley   } else PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, residualctx));
17359566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, jacobianctx));
17369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex));
17373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17389f520fc2SToby Isaac }
17399f520fc2SToby Isaac 
1740f017bcb6SMatthew G. Knepley /*@C
1741f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1742f017bcb6SMatthew G. Knepley 
1743f017bcb6SMatthew G. Knepley   Input Parameters:
1744f6dfbefdSBarry Smith + snes - the `SNES` object
1745f6dfbefdSBarry Smith . dm   - the `DM`
1746f2cacb80SMatthew G. Knepley . t    - the time
1747f6dfbefdSBarry Smith . u    - a `DM` vector
1748f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1749f017bcb6SMatthew G. Knepley 
17502fe279fdSBarry Smith   Output Parameter:
17512fe279fdSBarry Smith . error - An array which holds the discretization error in each field, or `NULL`
17522fe279fdSBarry Smith 
17532fe279fdSBarry Smith   Level: developer
1754f3c94560SMatthew G. Knepley 
1755f6dfbefdSBarry Smith   Note:
1756f6dfbefdSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
17577f96f943SMatthew G. Knepley 
1758e4094ef1SJacob Faibussowitsch .seealso: `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`
1759f017bcb6SMatthew G. Knepley @*/
1760d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
1761d71ae5a4SJacob Faibussowitsch {
1762f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1763f017bcb6SMatthew G. Knepley   void     **ectxs;
1764f3c94560SMatthew G. Knepley   PetscReal *err;
17657f96f943SMatthew G. Knepley   MPI_Comm   comm;
17667f96f943SMatthew G. Knepley   PetscInt   Nf, f;
17671878804aSMatthew G. Knepley 
17681878804aSMatthew G. Knepley   PetscFunctionBegin;
1769f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1770f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1771064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
17724f572ea9SToby Isaac   if (error) PetscAssertPointer(error, 6);
17737f96f943SMatthew G. Knepley 
17749566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
17759566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
17767f96f943SMatthew G. Knepley 
17779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
17789566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
17799566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
17807f96f943SMatthew G. Knepley   {
17817f96f943SMatthew G. Knepley     PetscInt Nds, s;
17827f96f943SMatthew G. Knepley 
17839566063dSJacob Faibussowitsch     PetscCall(DMGetNumDS(dm, &Nds));
1784083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
17857f96f943SMatthew G. Knepley       PetscDS         ds;
1786083401c6SMatthew G. Knepley       DMLabel         label;
1787083401c6SMatthew G. Knepley       IS              fieldIS;
17887f96f943SMatthew G. Knepley       const PetscInt *fields;
17897f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1790083401c6SMatthew G. Knepley 
179107218a29SMatthew G. Knepley       PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
17929566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(ds, &dsNf));
17939566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
1794083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1795083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
17969566063dSJacob Faibussowitsch         PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
1797083401c6SMatthew G. Knepley       }
17989566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
1799083401c6SMatthew G. Knepley     }
1800083401c6SMatthew G. Knepley   }
1801f017bcb6SMatthew G. Knepley   if (Nf > 1) {
18029566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err));
1803f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1804ad540459SPierre 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);
1805b878532fSMatthew G. Knepley     } else if (error) {
1806b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
18071878804aSMatthew G. Knepley     } else {
18089566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: ["));
1809f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
18109566063dSJacob Faibussowitsch         if (f) PetscCall(PetscPrintf(comm, ", "));
18119566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(comm, "%g", (double)err[f]));
18121878804aSMatthew G. Knepley       }
18139566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "]\n"));
1814f017bcb6SMatthew G. Knepley     }
1815f017bcb6SMatthew G. Knepley   } else {
18169566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]));
1817f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
181808401ef6SPierre Jolivet       PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol);
1819b878532fSMatthew G. Knepley     } else if (error) {
1820b878532fSMatthew G. Knepley       error[0] = err[0];
1821f017bcb6SMatthew G. Knepley     } else {
18229566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]));
1823f017bcb6SMatthew G. Knepley     }
1824f017bcb6SMatthew G. Knepley   }
18259566063dSJacob Faibussowitsch   PetscCall(PetscFree3(exacts, ectxs, err));
18263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1827f017bcb6SMatthew G. Knepley }
1828f017bcb6SMatthew G. Knepley 
1829f017bcb6SMatthew G. Knepley /*@C
1830f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1831f017bcb6SMatthew G. Knepley 
1832f017bcb6SMatthew G. Knepley   Input Parameters:
1833f6dfbefdSBarry Smith + snes - the `SNES` object
1834f6dfbefdSBarry Smith . dm   - the `DM`
1835f6dfbefdSBarry Smith . u    - a `DM` vector
1836f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1837f017bcb6SMatthew G. Knepley 
1838f6dfbefdSBarry Smith   Output Parameter:
18392fe279fdSBarry Smith . residual - The residual norm of the exact solution, or `NULL`
1840f3c94560SMatthew G. Knepley 
1841f017bcb6SMatthew G. Knepley   Level: developer
1842f017bcb6SMatthew G. Knepley 
1843db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
1844f017bcb6SMatthew G. Knepley @*/
1845d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1846d71ae5a4SJacob Faibussowitsch {
1847f017bcb6SMatthew G. Knepley   MPI_Comm  comm;
1848f017bcb6SMatthew G. Knepley   Vec       r;
1849f017bcb6SMatthew G. Knepley   PetscReal res;
1850f017bcb6SMatthew G. Knepley 
1851f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1852f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1853f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1854f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
18554f572ea9SToby Isaac   if (residual) PetscAssertPointer(residual, 5);
18569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
18579566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
18589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
18599566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, u, r));
18609566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
1861f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
186208401ef6SPierre Jolivet     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
1863b878532fSMatthew G. Knepley   } else if (residual) {
1864b878532fSMatthew G. Knepley     *residual = res;
1865f017bcb6SMatthew G. Knepley   } else {
18669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
186793ec0da9SPierre Jolivet     PetscCall(VecFilter(r, 1.0e-10));
18689566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
18699566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
18709566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
1871f017bcb6SMatthew G. Knepley   }
18729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
18733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1874f017bcb6SMatthew G. Knepley }
1875f017bcb6SMatthew G. Knepley 
1876f017bcb6SMatthew G. Knepley /*@C
1877f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1878f017bcb6SMatthew G. Knepley 
1879f017bcb6SMatthew G. Knepley   Input Parameters:
1880f6dfbefdSBarry Smith + snes - the `SNES` object
1881f6dfbefdSBarry Smith . dm   - the `DM`
1882f6dfbefdSBarry Smith . u    - a `DM` vector
1883f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1884f017bcb6SMatthew G. Knepley 
1885f3c94560SMatthew G. Knepley   Output Parameters:
18862fe279fdSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL`
18872fe279fdSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL`
1888f3c94560SMatthew G. Knepley 
1889f017bcb6SMatthew G. Knepley   Level: developer
1890f017bcb6SMatthew G. Knepley 
1891db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
1892f017bcb6SMatthew G. Knepley @*/
1893d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1894d71ae5a4SJacob Faibussowitsch {
1895f017bcb6SMatthew G. Knepley   MPI_Comm     comm;
1896f017bcb6SMatthew G. Knepley   PetscDS      ds;
1897f017bcb6SMatthew G. Knepley   Mat          J, M;
1898f017bcb6SMatthew G. Knepley   MatNullSpace nullspace;
1899f3c94560SMatthew G. Knepley   PetscReal    slope, intercept;
19007f96f943SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
1901f017bcb6SMatthew G. Knepley 
1902f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1903f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1904f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1905f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
19064f572ea9SToby Isaac   if (isLinear) PetscAssertPointer(isLinear, 5);
19074f572ea9SToby Isaac   if (convRate) PetscAssertPointer(convRate, 6);
19089566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
19099566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
1910f017bcb6SMatthew G. Knepley   /* Create and view matrices */
19119566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
19129566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
19139566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
19149566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
1915282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
19169566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
19179566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, M));
19189566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
19199566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
19209566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
19219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
1922282e7bb4SMatthew G. Knepley   } else {
19239566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, J));
1924282e7bb4SMatthew G. Knepley   }
19259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
19269566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
19279566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
1928f017bcb6SMatthew G. Knepley   /* Check nullspace */
19299566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
1930f017bcb6SMatthew G. Knepley   if (nullspace) {
19311878804aSMatthew G. Knepley     PetscBool isNull;
19329566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
193328b400f6SJacob Faibussowitsch     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
19341878804aSMatthew G. Knepley   }
1935f017bcb6SMatthew G. Knepley   /* Taylor test */
1936f017bcb6SMatthew G. Knepley   {
1937f017bcb6SMatthew G. Knepley     PetscRandom rand;
1938f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1939f3c94560SMatthew G. Knepley     PetscReal   h;
1940f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1941f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1942f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1943f017bcb6SMatthew G. Knepley 
1944f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
19459566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
19469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
19479566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
19489566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
19499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
19509566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
1951f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
19529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
19539566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, u, r));
1954f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
19559371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
19569371c9d4SSatish Balay       ;
19579566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
19589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
19599566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
1960f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
19619566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
1962f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
19639566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, uhat, rhat));
19649566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
19659566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
1966f017bcb6SMatthew G. Knepley 
196738d69901SMatthew G. Knepley       es[Nv] = errors[Nv] == 0 ? -16.0 : PetscLog10Real(errors[Nv]);
1968f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1969f017bcb6SMatthew G. Knepley     }
19709566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
19719566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
19729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
19739566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
19749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
1975f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1976f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1977f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
1978f017bcb6SMatthew G. Knepley     }
1979f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
19809566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
19819566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
1982f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1983f017bcb6SMatthew G. Knepley     if (tol >= 0) {
19840b121fc5SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
1985b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1986b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1987b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1988f017bcb6SMatthew G. Knepley     } else {
19899566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
19909566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
1991f017bcb6SMatthew G. Knepley     }
1992f017bcb6SMatthew G. Knepley   }
19939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
19943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19951878804aSMatthew G. Knepley }
19961878804aSMatthew G. Knepley 
1997d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1998d71ae5a4SJacob Faibussowitsch {
1999f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
20009566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL));
20019566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL));
20029566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL));
20033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2004f017bcb6SMatthew G. Knepley }
2005f017bcb6SMatthew G. Knepley 
2006bee9a294SMatthew G. Knepley /*@C
2007bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
2008bee9a294SMatthew G. Knepley 
2009bee9a294SMatthew G. Knepley   Input Parameters:
2010f6dfbefdSBarry Smith + snes - the `SNES` object
2011f6dfbefdSBarry Smith - u    - representative `SNES` vector
20127f96f943SMatthew G. Knepley 
20132fe279fdSBarry Smith   Level: developer
20142fe279fdSBarry Smith 
2015f6dfbefdSBarry Smith   Note:
2016f6dfbefdSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
2017bee9a294SMatthew G. Knepley 
20182fe279fdSBarry Smith .seealso: `SNES`, `DM`
2019bee9a294SMatthew G. Knepley @*/
2020d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
2021d71ae5a4SJacob Faibussowitsch {
20221878804aSMatthew G. Knepley   DM        dm;
20231878804aSMatthew G. Knepley   Vec       sol;
20241878804aSMatthew G. Knepley   PetscBool check;
20251878804aSMatthew G. Knepley 
20261878804aSMatthew G. Knepley   PetscFunctionBegin;
20279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check));
20283ba16761SJacob Faibussowitsch   if (!check) PetscFunctionReturn(PETSC_SUCCESS);
20299566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
20309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
20319566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, sol));
20329566063dSJacob Faibussowitsch   PetscCall(DMSNESCheck_Internal(snes, dm, sol));
20339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
20343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2035552f7358SJed Brown }
2036