xref: /petsc/src/snes/utils/dmplexsnes.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
79371c9d4SSatish Balay 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[]) {
8649ef022SMatthew Knepley   p[0] = u[uOff[1]];
9649ef022SMatthew Knepley }
10649ef022SMatthew Knepley 
11649ef022SMatthew Knepley /*
12649ef022SMatthew Knepley   SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately.
13649ef022SMatthew Knepley 
14*f6dfbefdSBarry Smith   Collective on snes
15649ef022SMatthew Knepley 
16649ef022SMatthew Knepley   Input Parameters:
17649ef022SMatthew Knepley + snes      - The SNES
18649ef022SMatthew Knepley . pfield    - The field number for pressure
19649ef022SMatthew Knepley . nullspace - The pressure nullspace
20649ef022SMatthew Knepley . u         - The solution vector
21649ef022SMatthew Knepley - ctx       - An optional user context
22649ef022SMatthew Knepley 
23649ef022SMatthew Knepley   Output Parameter:
24649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
25649ef022SMatthew Knepley 
26649ef022SMatthew Knepley   Notes:
27649ef022SMatthew 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.
28649ef022SMatthew Knepley 
29649ef022SMatthew Knepley   Level: developer
30649ef022SMatthew Knepley 
31db781477SPatrick Sanan .seealso: `SNESConvergedCorrectPressure()`
32649ef022SMatthew Knepley */
339371c9d4SSatish Balay static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx) {
34649ef022SMatthew Knepley   DM          dm;
35649ef022SMatthew Knepley   PetscDS     ds;
36649ef022SMatthew Knepley   const Vec  *nullvecs;
37649ef022SMatthew Knepley   PetscScalar pintd, *intc, *intn;
38649ef022SMatthew Knepley   MPI_Comm    comm;
39649ef022SMatthew Knepley   PetscInt    Nf, Nv;
40649ef022SMatthew Knepley 
41649ef022SMatthew Knepley   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
439566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4428b400f6SJacob Faibussowitsch   PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
4528b400f6SJacob Faibussowitsch   PetscCheck(nullspace, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
469566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
479566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, pfield, pressure_Private));
489566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs));
4963a3b9bcSJacob Faibussowitsch   PetscCheck(Nv == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %" PetscInt_FMT, Nv);
509566063dSJacob Faibussowitsch   PetscCall(VecDot(nullvecs[0], u, &pintd));
5108401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double)PetscRealPart(pintd));
529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &intc, Nf, &intn));
549566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx));
559566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
569566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u, -intc[pfield] / intn[pfield], nullvecs[0]));
57649ef022SMatthew Knepley #if defined(PETSC_USE_DEBUG)
589566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
5908401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(intc[pfield]) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double)PetscRealPart(intc[pfield]));
60649ef022SMatthew Knepley #endif
619566063dSJacob Faibussowitsch   PetscCall(PetscFree2(intc, intn));
62649ef022SMatthew Knepley   PetscFunctionReturn(0);
63649ef022SMatthew Knepley }
64649ef022SMatthew Knepley 
65649ef022SMatthew Knepley /*@C
66*f6dfbefdSBarry Smith    SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero.
67*f6dfbefdSBarry Smith    This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics `SNESConvergedDefault()`.
68649ef022SMatthew Knepley 
69*f6dfbefdSBarry Smith    Logically Collective on snes
70649ef022SMatthew Knepley 
71649ef022SMatthew Knepley    Input Parameters:
72*f6dfbefdSBarry Smith +  snes - the `SNES` context
73649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
74649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
75649ef022SMatthew Knepley .  snorm - 2-norm of current step
76649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
77649ef022SMatthew Knepley -  ctx   - Optional user context
78649ef022SMatthew Knepley 
79649ef022SMatthew Knepley    Output Parameter:
80*f6dfbefdSBarry Smith .  reason  - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FNORM_NAN`
81649ef022SMatthew Knepley 
82649ef022SMatthew Knepley    Notes:
83*f6dfbefdSBarry 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`
84*f6dfbefdSBarry Smith    must be created with discretizations of those fields. We currently assume that the pressure field has index 1.
85*f6dfbefdSBarry Smith    The pressure field must have a nullspace, likely created using the `DMSetNullSpaceConstructor()` interface.
86*f6dfbefdSBarry 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.
87f362779dSJed Brown 
88*f6dfbefdSBarry Smith    Options Database Key:
89*f6dfbefdSBarry Smith .  -snes_convergence_test correct_pressure - see `SNESSetFromOptions()`
90649ef022SMatthew Knepley 
91649ef022SMatthew Knepley    Level: advanced
92649ef022SMatthew Knepley 
93*f6dfbefdSBarry Smith .seealso: `SNES`, `DM`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()`, `DMSetNullSpaceConstructor()`
94649ef022SMatthew Knepley @*/
959371c9d4SSatish Balay PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx) {
96649ef022SMatthew Knepley   PetscBool monitorIntegral = PETSC_FALSE;
97649ef022SMatthew Knepley 
98649ef022SMatthew Knepley   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx));
100649ef022SMatthew Knepley   if (monitorIntegral) {
101649ef022SMatthew Knepley     Mat          J;
102649ef022SMatthew Knepley     Vec          u;
103649ef022SMatthew Knepley     MatNullSpace nullspace;
104649ef022SMatthew Knepley     const Vec   *nullvecs;
105649ef022SMatthew Knepley     PetscScalar  pintd;
106649ef022SMatthew Knepley 
1079566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1089566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1099566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1109566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs));
1119566063dSJacob Faibussowitsch     PetscCall(VecDot(nullvecs[0], u, &pintd));
1129566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double)PetscRealPart(pintd)));
113649ef022SMatthew Knepley   }
114649ef022SMatthew Knepley   if (*reason > 0) {
115649ef022SMatthew Knepley     Mat          J;
116649ef022SMatthew Knepley     Vec          u;
117649ef022SMatthew Knepley     MatNullSpace nullspace;
118649ef022SMatthew Knepley     PetscInt     pfield = 1;
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(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx));
124649ef022SMatthew Knepley   }
125649ef022SMatthew Knepley   PetscFunctionReturn(0);
126649ef022SMatthew Knepley }
127649ef022SMatthew Knepley 
12824cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
12924cdb843SMatthew G. Knepley 
1309371c9d4SSatish Balay static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) {
1316da023fcSToby Isaac   PetscBool isPlex;
1326da023fcSToby Isaac 
1336da023fcSToby Isaac   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
1356da023fcSToby Isaac   if (isPlex) {
1366da023fcSToby Isaac     *plex = dm;
1379566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
138f7148743SMatthew G. Knepley   } else {
1399566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
140f7148743SMatthew G. Knepley     if (!*plex) {
1419566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
1429566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
1436da023fcSToby Isaac       if (copy) {
1449566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
1459566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
1466da023fcSToby Isaac       }
147f7148743SMatthew G. Knepley     } else {
1489566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*plex));
149f7148743SMatthew G. Knepley     }
1506da023fcSToby Isaac   }
1516da023fcSToby Isaac   PetscFunctionReturn(0);
1526da023fcSToby Isaac }
1536da023fcSToby Isaac 
1544267b1a3SMatthew G. Knepley /*@C
155*f6dfbefdSBarry Smith   DMInterpolationCreate - Creates a `DMInterpolationInfo` context
1564267b1a3SMatthew G. Knepley 
157d083f849SBarry Smith   Collective
1584267b1a3SMatthew G. Knepley 
1594267b1a3SMatthew G. Knepley   Input Parameter:
1604267b1a3SMatthew G. Knepley . comm - the communicator
1614267b1a3SMatthew G. Knepley 
1624267b1a3SMatthew G. Knepley   Output Parameter:
1634267b1a3SMatthew G. Knepley . ctx - the context
1644267b1a3SMatthew G. Knepley 
1654267b1a3SMatthew G. Knepley   Level: beginner
1664267b1a3SMatthew G. Knepley 
167*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()`
1684267b1a3SMatthew G. Knepley @*/
1699371c9d4SSatish Balay PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) {
170552f7358SJed Brown   PetscFunctionBegin;
171552f7358SJed Brown   PetscValidPointer(ctx, 2);
1729566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
1731aa26658SKarl Rupp 
174552f7358SJed Brown   (*ctx)->comm   = comm;
175552f7358SJed Brown   (*ctx)->dim    = -1;
176552f7358SJed Brown   (*ctx)->nInput = 0;
1770298fd71SBarry Smith   (*ctx)->points = NULL;
1780298fd71SBarry Smith   (*ctx)->cells  = NULL;
179552f7358SJed Brown   (*ctx)->n      = -1;
1800298fd71SBarry Smith   (*ctx)->coords = NULL;
181552f7358SJed Brown   PetscFunctionReturn(0);
182552f7358SJed Brown }
183552f7358SJed Brown 
1844267b1a3SMatthew G. Knepley /*@C
1854267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
1864267b1a3SMatthew G. Knepley 
1874267b1a3SMatthew G. Knepley   Not collective
1884267b1a3SMatthew G. Knepley 
1894267b1a3SMatthew G. Knepley   Input Parameters:
1904267b1a3SMatthew G. Knepley + ctx - the context
1914267b1a3SMatthew G. Knepley - dim - the spatial dimension
1924267b1a3SMatthew G. Knepley 
1934267b1a3SMatthew G. Knepley   Level: intermediate
1944267b1a3SMatthew G. Knepley 
195*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
1964267b1a3SMatthew G. Knepley @*/
1979371c9d4SSatish Balay PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) {
198552f7358SJed Brown   PetscFunctionBegin;
19963a3b9bcSJacob Faibussowitsch   PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim);
200552f7358SJed Brown   ctx->dim = dim;
201552f7358SJed Brown   PetscFunctionReturn(0);
202552f7358SJed Brown }
203552f7358SJed Brown 
2044267b1a3SMatthew G. Knepley /*@C
2054267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2064267b1a3SMatthew G. Knepley 
2074267b1a3SMatthew G. Knepley   Not collective
2084267b1a3SMatthew G. Knepley 
2094267b1a3SMatthew G. Knepley   Input Parameter:
2104267b1a3SMatthew G. Knepley . ctx - the context
2114267b1a3SMatthew G. Knepley 
2124267b1a3SMatthew G. Knepley   Output Parameter:
2134267b1a3SMatthew G. Knepley . dim - the spatial dimension
2144267b1a3SMatthew G. Knepley 
2154267b1a3SMatthew G. Knepley   Level: intermediate
2164267b1a3SMatthew G. Knepley 
217*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2184267b1a3SMatthew G. Knepley @*/
2199371c9d4SSatish Balay PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) {
220552f7358SJed Brown   PetscFunctionBegin;
221552f7358SJed Brown   PetscValidIntPointer(dim, 2);
222552f7358SJed Brown   *dim = ctx->dim;
223552f7358SJed Brown   PetscFunctionReturn(0);
224552f7358SJed Brown }
225552f7358SJed Brown 
2264267b1a3SMatthew G. Knepley /*@C
2274267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2284267b1a3SMatthew G. Knepley 
2294267b1a3SMatthew G. Knepley   Not collective
2304267b1a3SMatthew G. Knepley 
2314267b1a3SMatthew G. Knepley   Input Parameters:
2324267b1a3SMatthew G. Knepley + ctx - the context
2334267b1a3SMatthew G. Knepley - dof - the number of fields
2344267b1a3SMatthew G. Knepley 
2354267b1a3SMatthew G. Knepley   Level: intermediate
2364267b1a3SMatthew G. Knepley 
237*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2384267b1a3SMatthew G. Knepley @*/
2399371c9d4SSatish Balay PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) {
240552f7358SJed Brown   PetscFunctionBegin;
24163a3b9bcSJacob Faibussowitsch   PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof);
242552f7358SJed Brown   ctx->dof = dof;
243552f7358SJed Brown   PetscFunctionReturn(0);
244552f7358SJed Brown }
245552f7358SJed Brown 
2464267b1a3SMatthew G. Knepley /*@C
2474267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2484267b1a3SMatthew G. Knepley 
2494267b1a3SMatthew G. Knepley   Not collective
2504267b1a3SMatthew G. Knepley 
2514267b1a3SMatthew G. Knepley   Input Parameter:
2524267b1a3SMatthew G. Knepley . ctx - the context
2534267b1a3SMatthew G. Knepley 
2544267b1a3SMatthew G. Knepley   Output Parameter:
2554267b1a3SMatthew G. Knepley . dof - the number of fields
2564267b1a3SMatthew G. Knepley 
2574267b1a3SMatthew G. Knepley   Level: intermediate
2584267b1a3SMatthew G. Knepley 
259*f6dfbefdSBarry Smith .seealso: DMInterpolationInfo, `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2604267b1a3SMatthew G. Knepley @*/
2619371c9d4SSatish Balay PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) {
262552f7358SJed Brown   PetscFunctionBegin;
263552f7358SJed Brown   PetscValidIntPointer(dof, 2);
264552f7358SJed Brown   *dof = ctx->dof;
265552f7358SJed Brown   PetscFunctionReturn(0);
266552f7358SJed Brown }
267552f7358SJed Brown 
2684267b1a3SMatthew G. Knepley /*@C
2694267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2704267b1a3SMatthew G. Knepley 
2714267b1a3SMatthew G. Knepley   Not collective
2724267b1a3SMatthew G. Knepley 
2734267b1a3SMatthew G. Knepley   Input Parameters:
2744267b1a3SMatthew G. Knepley + ctx    - the context
2754267b1a3SMatthew G. Knepley . n      - the number of points
2764267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2774267b1a3SMatthew G. Knepley 
278*f6dfbefdSBarry Smith   Note:
279*f6dfbefdSBarry Smith   The coordinate information is copied.
2804267b1a3SMatthew G. Knepley 
2814267b1a3SMatthew G. Knepley   Level: intermediate
2824267b1a3SMatthew G. Knepley 
283*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()`
2844267b1a3SMatthew G. Knepley @*/
2859371c9d4SSatish Balay PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) {
286552f7358SJed Brown   PetscFunctionBegin;
28708401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
28828b400f6SJacob Faibussowitsch   PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
289552f7358SJed Brown   ctx->nInput = n;
2901aa26658SKarl Rupp 
2919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points));
2929566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim));
293552f7358SJed Brown   PetscFunctionReturn(0);
294552f7358SJed Brown }
295552f7358SJed Brown 
2964267b1a3SMatthew G. Knepley /*@C
29752aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
2984267b1a3SMatthew G. Knepley 
2994267b1a3SMatthew G. Knepley   Collective on ctx
3004267b1a3SMatthew G. Knepley 
3014267b1a3SMatthew G. Knepley   Input Parameters:
3024267b1a3SMatthew G. Knepley + ctx - the context
303*f6dfbefdSBarry Smith . dm  - the `DM` for the function space used for interpolation
304*f6dfbefdSBarry Smith . redundantPoints - If `PETSC_TRUE`, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
305*f6dfbefdSBarry Smith - ignoreOutsideDomain - If `PETSC_TRUE`, ignore points outside the domain, otherwise return an error
3064267b1a3SMatthew G. Knepley 
3074267b1a3SMatthew G. Knepley   Level: intermediate
3084267b1a3SMatthew G. Knepley 
309*f6dfbefdSBarry Smith .seealso: DMInterpolationInfo, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
3104267b1a3SMatthew G. Knepley @*/
3119371c9d4SSatish Balay PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain) {
312552f7358SJed Brown   MPI_Comm           comm = ctx->comm;
313552f7358SJed Brown   PetscScalar       *a;
314552f7358SJed Brown   PetscInt           p, q, i;
315552f7358SJed Brown   PetscMPIInt        rank, size;
316552f7358SJed Brown   Vec                pointVec;
3173a93e3b7SToby Isaac   PetscSF            cellSF;
318552f7358SJed Brown   PetscLayout        layout;
319552f7358SJed Brown   PetscReal         *globalPoints;
320cb313848SJed Brown   PetscScalar       *globalPointsScalar;
321552f7358SJed Brown   const PetscInt    *ranges;
322552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3233a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3243a93e3b7SToby Isaac   const PetscInt    *foundPoints;
325552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3263a93e3b7SToby Isaac   PetscInt           n, N, numFound;
327552f7358SJed Brown 
32819436ca2SJed Brown   PetscFunctionBegin;
329064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
3309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
33208401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
33319436ca2SJed Brown   /* Locate points */
33419436ca2SJed Brown   n = ctx->nInput;
335552f7358SJed Brown   if (!redundantPoints) {
3369566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &layout));
3379566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(layout, 1));
3389566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(layout, n));
3399566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(layout));
3409566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(layout, &N));
341552f7358SJed Brown     /* Communicate all points to all processes */
3429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs));
3439566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(layout, &ranges));
344552f7358SJed Brown     for (p = 0; p < size; ++p) {
345552f7358SJed Brown       counts[p] = (ranges[p + 1] - ranges[p]) * ctx->dim;
346552f7358SJed Brown       displs[p] = ranges[p] * ctx->dim;
347552f7358SJed Brown     }
3489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(ctx->points, n * ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm));
349552f7358SJed Brown   } else {
350552f7358SJed Brown     N            = n;
351552f7358SJed Brown     globalPoints = ctx->points;
35238ea73c8SJed Brown     counts = displs = NULL;
35338ea73c8SJed Brown     layout          = NULL;
354552f7358SJed Brown   }
355552f7358SJed Brown #if 0
3569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs));
35719436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
358552f7358SJed Brown #else
359cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
3609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar));
36146b3086cSToby Isaac   for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
362cb313848SJed Brown #else
363cb313848SJed Brown   globalPointsScalar = globalPoints;
364cb313848SJed Brown #endif
3659566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec));
3669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs));
367ad540459SPierre Jolivet   for (p = 0; p < N; ++p) foundProcs[p] = size;
3683a93e3b7SToby Isaac   cellSF = NULL;
3699566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF));
3709566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells));
371552f7358SJed Brown #endif
3723a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3733a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
374552f7358SJed Brown   }
375552f7358SJed Brown   /* Let the lowest rank process own each point */
3761c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm));
377552f7358SJed Brown   ctx->n = 0;
378552f7358SJed Brown   for (p = 0; p < N; ++p) {
37952aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
3809371c9d4SSatish 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),
3819371c9d4SSatish Balay                  (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0.0));
382f7d195e4SLawrence Mitchell       if (rank == 0) ++ctx->n;
38352aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
384552f7358SJed Brown   }
385552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
3869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->n, &ctx->cells));
3879566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &ctx->coords));
3889566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE));
3899566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(ctx->coords, ctx->dim));
3909566063dSJacob Faibussowitsch   PetscCall(VecSetType(ctx->coords, VECSTANDARD));
3919566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->coords, &a));
392552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
393552f7358SJed Brown     if (globalProcs[p] == rank) {
394552f7358SJed Brown       PetscInt d;
395552f7358SJed Brown 
3961aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d];
397f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
398f317b747SMatthew G. Knepley       ++q;
399552f7358SJed Brown     }
400dd400576SPatrick Sanan     if (globalProcs[p] == size && rank == 0) {
40152aa1562SMatthew G. Knepley       PetscInt d;
40252aa1562SMatthew G. Knepley 
40352aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
40452aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
40552aa1562SMatthew G. Knepley       ++q;
40652aa1562SMatthew G. Knepley     }
407552f7358SJed Brown   }
4089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->coords, &a));
409552f7358SJed Brown #if 0
4109566063dSJacob Faibussowitsch   PetscCall(PetscFree3(foundCells,foundProcs,globalProcs));
411552f7358SJed Brown #else
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(foundProcs, globalProcs));
4139566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
4149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pointVec));
415552f7358SJed Brown #endif
4169566063dSJacob Faibussowitsch   if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar));
4179566063dSJacob Faibussowitsch   if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs));
4189566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
419552f7358SJed Brown   PetscFunctionReturn(0);
420552f7358SJed Brown }
421552f7358SJed Brown 
4224267b1a3SMatthew G. Knepley /*@C
423*f6dfbefdSBarry Smith   DMInterpolationGetCoordinates - Gets a `Vec` with the coordinates of each interpolation point
4244267b1a3SMatthew G. Knepley 
4254267b1a3SMatthew G. Knepley   Collective on ctx
4264267b1a3SMatthew G. Knepley 
4274267b1a3SMatthew G. Knepley   Input Parameter:
4284267b1a3SMatthew G. Knepley . ctx - the context
4294267b1a3SMatthew G. Knepley 
4304267b1a3SMatthew G. Knepley   Output Parameter:
4314267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4324267b1a3SMatthew G. Knepley 
433*f6dfbefdSBarry Smith   Note:
434*f6dfbefdSBarry Smith   The local vector entries correspond to interpolation points lying on this process, according to the associated `DM`.
435*f6dfbefdSBarry Smith   This is a borrowed vector that the user should not destroy.
4364267b1a3SMatthew G. Knepley 
4374267b1a3SMatthew G. Knepley   Level: intermediate
4384267b1a3SMatthew G. Knepley 
439*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4404267b1a3SMatthew G. Knepley @*/
4419371c9d4SSatish Balay PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) {
442552f7358SJed Brown   PetscFunctionBegin;
443552f7358SJed Brown   PetscValidPointer(coordinates, 2);
44428b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
445552f7358SJed Brown   *coordinates = ctx->coords;
446552f7358SJed Brown   PetscFunctionReturn(0);
447552f7358SJed Brown }
448552f7358SJed Brown 
4494267b1a3SMatthew G. Knepley /*@C
450*f6dfbefdSBarry Smith   DMInterpolationGetVector - Gets a `Vec` which can hold all the interpolated field values
4514267b1a3SMatthew G. Knepley 
4524267b1a3SMatthew G. Knepley   Collective on ctx
4534267b1a3SMatthew G. Knepley 
4544267b1a3SMatthew G. Knepley   Input Parameter:
4554267b1a3SMatthew G. Knepley . ctx - the context
4564267b1a3SMatthew G. Knepley 
4574267b1a3SMatthew G. Knepley   Output Parameter:
4584267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4594267b1a3SMatthew G. Knepley 
460*f6dfbefdSBarry Smith   Note:
461*f6dfbefdSBarry Smith   This vector should be returned using `DMInterpolationRestoreVector()`.
4624267b1a3SMatthew G. Knepley 
4634267b1a3SMatthew G. Knepley   Level: intermediate
4644267b1a3SMatthew G. Knepley 
465*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4664267b1a3SMatthew G. Knepley @*/
4679371c9d4SSatish Balay PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) {
468552f7358SJed Brown   PetscFunctionBegin;
469552f7358SJed Brown   PetscValidPointer(v, 2);
47028b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4719566063dSJacob Faibussowitsch   PetscCall(VecCreate(ctx->comm, v));
4729566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE));
4739566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*v, ctx->dof));
4749566063dSJacob Faibussowitsch   PetscCall(VecSetType(*v, VECSTANDARD));
475552f7358SJed Brown   PetscFunctionReturn(0);
476552f7358SJed Brown }
477552f7358SJed Brown 
4784267b1a3SMatthew G. Knepley /*@C
479*f6dfbefdSBarry Smith   DMInterpolationRestoreVector - Returns a `Vec` which can hold all the interpolated field values
4804267b1a3SMatthew G. Knepley 
4814267b1a3SMatthew G. Knepley   Collective on ctx
4824267b1a3SMatthew G. Knepley 
4834267b1a3SMatthew G. Knepley   Input Parameters:
4844267b1a3SMatthew G. Knepley + ctx - the context
4854267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
4864267b1a3SMatthew G. Knepley 
4874267b1a3SMatthew G. Knepley   Level: intermediate
4884267b1a3SMatthew G. Knepley 
489*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4904267b1a3SMatthew G. Knepley @*/
4919371c9d4SSatish Balay PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) {
492552f7358SJed Brown   PetscFunctionBegin;
493552f7358SJed Brown   PetscValidPointer(v, 2);
49428b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(v));
496552f7358SJed Brown   PetscFunctionReturn(0);
497552f7358SJed Brown }
498552f7358SJed Brown 
4999371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
500cfe5744fSMatthew G. Knepley   PetscReal          v0, J, invJ, detJ;
501cfe5744fSMatthew G. Knepley   const PetscInt     dof = ctx->dof;
502cfe5744fSMatthew G. Knepley   const PetscScalar *coords;
503cfe5744fSMatthew G. Knepley   PetscScalar       *a;
504cfe5744fSMatthew G. Knepley   PetscInt           p;
505cfe5744fSMatthew G. Knepley 
506cfe5744fSMatthew G. Knepley   PetscFunctionBegin;
5079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
509cfe5744fSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
510cfe5744fSMatthew G. Knepley     PetscInt     c = ctx->cells[p];
511cfe5744fSMatthew G. Knepley     PetscScalar *x = NULL;
512cfe5744fSMatthew G. Knepley     PetscReal    xir[1];
513cfe5744fSMatthew G. Knepley     PetscInt     xSize, comp;
514cfe5744fSMatthew G. Knepley 
5159566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ));
516cfe5744fSMatthew G. Knepley     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
517cfe5744fSMatthew G. Knepley     xir[0] = invJ * PetscRealPart(coords[p] - v0);
5189566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
519cfe5744fSMatthew G. Knepley     if (2 * dof == xSize) {
520cfe5744fSMatthew 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];
521cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
522cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
523cfe5744fSMatthew 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);
5249566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
525cfe5744fSMatthew G. Knepley   }
5269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
528cfe5744fSMatthew G. Knepley   PetscFunctionReturn(0);
529cfe5744fSMatthew G. Knepley }
530cfe5744fSMatthew G. Knepley 
5319371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
532552f7358SJed Brown   PetscReal         *v0, *J, *invJ, detJ;
53356044e6dSMatthew G. Knepley   const PetscScalar *coords;
53456044e6dSMatthew G. Knepley   PetscScalar       *a;
535552f7358SJed Brown   PetscInt           p;
536552f7358SJed Brown 
537552f7358SJed Brown   PetscFunctionBegin;
5389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5409566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
541552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
542552f7358SJed Brown     PetscInt     c = ctx->cells[p];
543a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
544552f7358SJed Brown     PetscReal    xi[4];
545552f7358SJed Brown     PetscInt     d, f, comp;
546552f7358SJed Brown 
5479566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
54863a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5499566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5501aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
5511aa26658SKarl Rupp 
552552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
553552f7358SJed Brown       xi[d] = 0.0;
5541aa26658SKarl Rupp       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]);
5551aa26658SKarl Rupp       for (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];
556552f7358SJed Brown     }
5579566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
558552f7358SJed Brown   }
5599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5619566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
562552f7358SJed Brown   PetscFunctionReturn(0);
563552f7358SJed Brown }
564552f7358SJed Brown 
5659371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
5667a1931ceSMatthew G. Knepley   PetscReal         *v0, *J, *invJ, detJ;
56756044e6dSMatthew G. Knepley   const PetscScalar *coords;
56856044e6dSMatthew G. Knepley   PetscScalar       *a;
5697a1931ceSMatthew G. Knepley   PetscInt           p;
5707a1931ceSMatthew G. Knepley 
5717a1931ceSMatthew G. Knepley   PetscFunctionBegin;
5729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5757a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5767a1931ceSMatthew G. Knepley     PetscInt       c        = ctx->cells[p];
5777a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5782584bbe8SMatthew G. Knepley     PetscScalar   *x        = NULL;
5797a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5807a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5817a1931ceSMatthew G. Knepley 
5829566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
58363a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5849566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5857a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
5867a1931ceSMatthew G. Knepley 
5877a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5887a1931ceSMatthew G. Knepley       xi[d] = 0.0;
5897a1931ceSMatthew G. Knepley       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]);
5907a1931ceSMatthew G. Knepley       for (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];
5917a1931ceSMatthew G. Knepley     }
5929566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
5937a1931ceSMatthew G. Knepley   }
5949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5969566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
5977a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
5987a1931ceSMatthew G. Knepley }
5997a1931ceSMatthew G. Knepley 
6009371c9d4SSatish Balay static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) {
601552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
602552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
603552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
604552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
605552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
606552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
607552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
608552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
609552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
610552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
611552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
612552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
613552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
614552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
615552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
61656044e6dSMatthew G. Knepley   const PetscScalar *ref;
61756044e6dSMatthew G. Knepley   PetscScalar       *real;
618552f7358SJed Brown 
619552f7358SJed Brown   PetscFunctionBegin;
6209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
6219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
622552f7358SJed Brown   {
623552f7358SJed Brown     const PetscScalar p0 = ref[0];
624552f7358SJed Brown     const PetscScalar p1 = ref[1];
625552f7358SJed Brown 
626552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
627552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
628552f7358SJed Brown   }
6299566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(28));
6309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
632552f7358SJed Brown   PetscFunctionReturn(0);
633552f7358SJed Brown }
634552f7358SJed Brown 
635af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
6369371c9d4SSatish Balay static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) {
637552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
638552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
639552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
640552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
641552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
642552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
643552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
644552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
645552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
646552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
647552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
64856044e6dSMatthew G. Knepley   const PetscScalar *ref;
649552f7358SJed Brown 
650552f7358SJed Brown   PetscFunctionBegin;
6519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
652552f7358SJed Brown   {
653552f7358SJed Brown     const PetscScalar x       = ref[0];
654552f7358SJed Brown     const PetscScalar y       = ref[1];
655552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
656da80777bSKarl Rupp     PetscScalar       values[4];
657da80777bSKarl Rupp 
6589371c9d4SSatish Balay     values[0] = (x1 - x0 + f_01 * y) * 0.5;
6599371c9d4SSatish Balay     values[1] = (x3 - x0 + f_01 * x) * 0.5;
6609371c9d4SSatish Balay     values[2] = (y1 - y0 + g_01 * y) * 0.5;
6619371c9d4SSatish Balay     values[3] = (y3 - y0 + g_01 * x) * 0.5;
6629566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES));
663552f7358SJed Brown   }
6649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(30));
6659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
668552f7358SJed Brown   PetscFunctionReturn(0);
669552f7358SJed Brown }
670552f7358SJed Brown 
6719371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
672fafc0619SMatthew G Knepley   DM                 dmCoord;
673de73a395SMatthew G. Knepley   PetscFE            fem = NULL;
674552f7358SJed Brown   SNES               snes;
675552f7358SJed Brown   KSP                ksp;
676552f7358SJed Brown   PC                 pc;
677552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
678552f7358SJed Brown   Mat                J;
679716009bfSMatthew G. Knepley   PetscTabulation    T = NULL;
68056044e6dSMatthew G. Knepley   const PetscScalar *coords;
68156044e6dSMatthew G. Knepley   PetscScalar       *a;
682716009bfSMatthew G. Knepley   PetscReal          xir[2] = {0., 0.};
683de73a395SMatthew G. Knepley   PetscInt           Nf, p;
6845509d985SMatthew G. Knepley   const PetscInt     dof = ctx->dof;
685552f7358SJed Brown 
686552f7358SJed Brown   PetscFunctionBegin;
6879566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
688716009bfSMatthew G. Knepley   if (Nf) {
689cfe5744fSMatthew G. Knepley     PetscObject  obj;
690cfe5744fSMatthew G. Knepley     PetscClassId id;
691cfe5744fSMatthew G. Knepley 
6929566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, 0, NULL, &obj));
6939566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
6949371c9d4SSatish Balay     if (id == PETSCFE_CLASSID) {
6959371c9d4SSatish Balay       fem = (PetscFE)obj;
6969371c9d4SSatish Balay       PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T));
6979371c9d4SSatish Balay     }
698716009bfSMatthew G. Knepley   }
6999566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
7009566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
7019566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
7029566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_"));
7039566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
7049566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 2, 2));
7059566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
7069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
7079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
7089566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
7099566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 2, 2, 2, 2));
7109566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
7119566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
7129566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL));
7139566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL));
7149566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7159566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
7169566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
7179566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
718552f7358SJed Brown 
7199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
7209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
721552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
722a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
723552f7358SJed Brown     PetscScalar *xi;
724552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
725552f7358SJed Brown 
726552f7358SJed Brown     /* Can make this do all points at once */
7279566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
72863a3b9bcSJacob Faibussowitsch     PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2);
7299566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
7309566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
7319566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
7329566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
733552f7358SJed Brown     xi[0] = coords[p * ctx->dim + 0];
734552f7358SJed Brown     xi[1] = coords[p * ctx->dim + 1];
7359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
7369566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
7379566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
738cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
739cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
740cfe5744fSMatthew G. Knepley     if (4 * dof == xSize) {
7419371c9d4SSatish Balay       for (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];
742cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
743cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
744cfe5744fSMatthew G. Knepley     } else {
7455509d985SMatthew G. Knepley       PetscInt d;
7461aa26658SKarl Rupp 
7472c71b3e2SJacob Faibussowitsch       PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE");
7489371c9d4SSatish Balay       xir[0] = 2.0 * xir[0] - 1.0;
7499371c9d4SSatish Balay       xir[1] = 2.0 * xir[1] - 1.0;
7509566063dSJacob Faibussowitsch       PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T));
7515509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7525509d985SMatthew G. Knepley         a[p * dof + comp] = 0.0;
753ad540459SPierre Jolivet         for (d = 0; d < xSize / dof; ++d) a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp];
7545509d985SMatthew G. Knepley       }
7555509d985SMatthew G. Knepley     }
7569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
7579566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
7589566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
759552f7358SJed Brown   }
7609566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&T));
7619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
7629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
763552f7358SJed Brown 
7649566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
7679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
7689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
769552f7358SJed Brown   PetscFunctionReturn(0);
770552f7358SJed Brown }
771552f7358SJed Brown 
7729371c9d4SSatish Balay static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) {
773552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
774552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
775552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
776552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
7777a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
7787a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
7797a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
780552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
781552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
782552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
7837a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
7847a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
7857a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
786552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
787552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
788552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
789552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
790552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
791552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
792552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
793552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
794552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
795552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
796552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
797552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
798552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
799552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
800552f7358SJed Brown   const PetscScalar  h_1      = z1 - z0;
801552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
802552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
803552f7358SJed Brown   const PetscScalar  h_3      = z3 - z0;
804552f7358SJed Brown   const PetscScalar  f_4      = x4 - x0;
805552f7358SJed Brown   const PetscScalar  g_4      = y4 - y0;
806552f7358SJed Brown   const PetscScalar  h_4      = z4 - z0;
807552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
808552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
809552f7358SJed Brown   const PetscScalar  h_01     = z2 - z1 - z3 + z0;
810552f7358SJed Brown   const PetscScalar  f_12     = x7 - x3 - x4 + x0;
811552f7358SJed Brown   const PetscScalar  g_12     = y7 - y3 - y4 + y0;
812552f7358SJed Brown   const PetscScalar  h_12     = z7 - z3 - z4 + z0;
813552f7358SJed Brown   const PetscScalar  f_02     = x5 - x1 - x4 + x0;
814552f7358SJed Brown   const PetscScalar  g_02     = y5 - y1 - y4 + y0;
815552f7358SJed Brown   const PetscScalar  h_02     = z5 - z1 - z4 + z0;
816552f7358SJed Brown   const PetscScalar  f_012    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
817552f7358SJed Brown   const PetscScalar  g_012    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
818552f7358SJed Brown   const PetscScalar  h_012    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
81956044e6dSMatthew G. Knepley   const PetscScalar *ref;
82056044e6dSMatthew G. Knepley   PetscScalar       *real;
821552f7358SJed Brown 
822552f7358SJed Brown   PetscFunctionBegin;
8239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
8249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
825552f7358SJed Brown   {
826552f7358SJed Brown     const PetscScalar p0 = ref[0];
827552f7358SJed Brown     const PetscScalar p1 = ref[1];
828552f7358SJed Brown     const PetscScalar p2 = ref[2];
829552f7358SJed Brown 
830552f7358SJed 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;
831552f7358SJed 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;
832552f7358SJed 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;
833552f7358SJed Brown   }
8349566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(114));
8359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
8369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
837552f7358SJed Brown   PetscFunctionReturn(0);
838552f7358SJed Brown }
839552f7358SJed Brown 
8409371c9d4SSatish Balay static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) {
841552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
842552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
843552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
844552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
8457a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
8467a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
8477a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
848552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
849552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
850552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
8517a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
8527a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
8537a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
854552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
855552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
856552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
857552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
858552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
859552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
860552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
861552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
862552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
863552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
864552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
865552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
866552f7358SJed Brown   const PetscScalar  f_xy     = x2 - x1 - x3 + x0;
867552f7358SJed Brown   const PetscScalar  g_xy     = y2 - y1 - y3 + y0;
868552f7358SJed Brown   const PetscScalar  h_xy     = z2 - z1 - z3 + z0;
869552f7358SJed Brown   const PetscScalar  f_yz     = x7 - x3 - x4 + x0;
870552f7358SJed Brown   const PetscScalar  g_yz     = y7 - y3 - y4 + y0;
871552f7358SJed Brown   const PetscScalar  h_yz     = z7 - z3 - z4 + z0;
872552f7358SJed Brown   const PetscScalar  f_xz     = x5 - x1 - x4 + x0;
873552f7358SJed Brown   const PetscScalar  g_xz     = y5 - y1 - y4 + y0;
874552f7358SJed Brown   const PetscScalar  h_xz     = z5 - z1 - z4 + z0;
875552f7358SJed Brown   const PetscScalar  f_xyz    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
876552f7358SJed Brown   const PetscScalar  g_xyz    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
877552f7358SJed Brown   const PetscScalar  h_xyz    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
87856044e6dSMatthew G. Knepley   const PetscScalar *ref;
879552f7358SJed Brown 
880552f7358SJed Brown   PetscFunctionBegin;
8819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
882552f7358SJed Brown   {
883552f7358SJed Brown     const PetscScalar x       = ref[0];
884552f7358SJed Brown     const PetscScalar y       = ref[1];
885552f7358SJed Brown     const PetscScalar z       = ref[2];
886552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
887da80777bSKarl Rupp     PetscScalar       values[9];
888da80777bSKarl Rupp 
889da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0;
890da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0;
891da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0;
892da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0;
893da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0;
894da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0;
895da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0;
896da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0;
897da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0;
8981aa26658SKarl Rupp 
8999566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES));
900552f7358SJed Brown   }
9019566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(152));
9029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
9039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
9049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
905552f7358SJed Brown   PetscFunctionReturn(0);
906552f7358SJed Brown }
907552f7358SJed Brown 
9089371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
909fafc0619SMatthew G Knepley   DM                 dmCoord;
910552f7358SJed Brown   SNES               snes;
911552f7358SJed Brown   KSP                ksp;
912552f7358SJed Brown   PC                 pc;
913552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
914552f7358SJed Brown   Mat                J;
91556044e6dSMatthew G. Knepley   const PetscScalar *coords;
91656044e6dSMatthew G. Knepley   PetscScalar       *a;
917552f7358SJed Brown   PetscInt           p;
918552f7358SJed Brown 
919552f7358SJed Brown   PetscFunctionBegin;
9209566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
9219566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
9229566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
9239566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_"));
9249566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
9259566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 3, 3));
9269566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
9279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
9289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
9299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
9309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 3, 3, 3, 3));
9319566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
9329566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
9339566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL));
9349566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL));
9359566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
9369566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
9379566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
9389566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
939552f7358SJed Brown 
9409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
9419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
942552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
943a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
944552f7358SJed Brown     PetscScalar *xi;
945cb313848SJed Brown     PetscReal    xir[3];
946552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
947552f7358SJed Brown 
948552f7358SJed Brown     /* Can make this do all points at once */
9499566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
950cfe5744fSMatthew G. Knepley     PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3);
9519566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
952cfe5744fSMatthew 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);
9539566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
9549566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
9559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
956552f7358SJed Brown     xi[0] = coords[p * ctx->dim + 0];
957552f7358SJed Brown     xi[1] = coords[p * ctx->dim + 1];
958552f7358SJed Brown     xi[2] = coords[p * ctx->dim + 2];
9599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
9609566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
9619566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
962cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
963cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
964cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
965cfe5744fSMatthew G. Knepley     if (8 * ctx->dof == xSize) {
966552f7358SJed Brown       for (comp = 0; comp < ctx->dof; ++comp) {
9679371c9d4SSatish 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]) +
9689371c9d4SSatish 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];
969552f7358SJed Brown       }
970cfe5744fSMatthew G. Knepley     } else {
971cfe5744fSMatthew G. Knepley       for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
972cfe5744fSMatthew G. Knepley     }
9739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
9749566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
9759566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
976552f7358SJed Brown   }
9779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
9789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
979552f7358SJed Brown 
9809566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
9819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
9829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
9839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
9849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
985552f7358SJed Brown   PetscFunctionReturn(0);
986552f7358SJed Brown }
987552f7358SJed Brown 
9884267b1a3SMatthew G. Knepley /*@C
9894267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9904267b1a3SMatthew G. Knepley 
991552f7358SJed Brown   Input Parameters:
992*f6dfbefdSBarry Smith + ctx - The `DMInterpolationInfo` context
993*f6dfbefdSBarry Smith . dm  - The `DM`
994552f7358SJed Brown - x   - The local vector containing the field to be interpolated
995552f7358SJed Brown 
996552f7358SJed Brown   Output Parameters:
997552f7358SJed Brown . v   - The vector containing the interpolated values
9984267b1a3SMatthew G. Knepley 
999*f6dfbefdSBarry Smith   Note:
1000*f6dfbefdSBarry Smith   A suitable v can be obtained using `DMInterpolationGetVector()`.
10014267b1a3SMatthew G. Knepley 
10024267b1a3SMatthew G. Knepley   Level: beginner
10034267b1a3SMatthew G. Knepley 
1004*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
10054267b1a3SMatthew G. Knepley @*/
10069371c9d4SSatish Balay PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) {
100766a0a883SMatthew G. Knepley   PetscDS   ds;
100866a0a883SMatthew G. Knepley   PetscInt  n, p, Nf, field;
100966a0a883SMatthew G. Knepley   PetscBool useDS = PETSC_FALSE;
1010552f7358SJed Brown 
1011552f7358SJed Brown   PetscFunctionBegin;
1012552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1013552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1014552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
10159566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
101663a3b9bcSJacob 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);
101766a0a883SMatthew G. Knepley   if (!n) PetscFunctionReturn(0);
10189566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1019680d18d5SMatthew G. Knepley   if (ds) {
102066a0a883SMatthew G. Knepley     useDS = PETSC_TRUE;
10219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &Nf));
1022680d18d5SMatthew G. Knepley     for (field = 0; field < Nf; ++field) {
102366a0a883SMatthew G. Knepley       PetscObject  obj;
1024680d18d5SMatthew G. Knepley       PetscClassId id;
1025680d18d5SMatthew G. Knepley 
10269566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, field, &obj));
10279566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
10289371c9d4SSatish Balay       if (id != PETSCFE_CLASSID) {
10299371c9d4SSatish Balay         useDS = PETSC_FALSE;
10309371c9d4SSatish Balay         break;
10319371c9d4SSatish Balay       }
103266a0a883SMatthew G. Knepley     }
103366a0a883SMatthew G. Knepley   }
103466a0a883SMatthew G. Knepley   if (useDS) {
103566a0a883SMatthew G. Knepley     const PetscScalar *coords;
103666a0a883SMatthew G. Knepley     PetscScalar       *interpolant;
103766a0a883SMatthew G. Knepley     PetscInt           cdim, d;
103866a0a883SMatthew G. Knepley 
10399566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm, &cdim));
10409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ctx->coords, &coords));
10419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(v, &interpolant));
1042680d18d5SMatthew G. Knepley     for (p = 0; p < ctx->n; ++p) {
104366a0a883SMatthew G. Knepley       PetscReal    pcoords[3], xi[3];
1044680d18d5SMatthew G. Knepley       PetscScalar *xa   = NULL;
104566a0a883SMatthew G. Knepley       PetscInt     coff = 0, foff = 0, clSize;
1046680d18d5SMatthew G. Knepley 
104752aa1562SMatthew G. Knepley       if (ctx->cells[p] < 0) continue;
1048680d18d5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]);
10499566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi));
10509566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
105166a0a883SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
105266a0a883SMatthew G. Knepley         PetscTabulation T;
105366a0a883SMatthew G. Knepley         PetscFE         fe;
105466a0a883SMatthew G. Knepley 
10559566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
10569566063dSJacob Faibussowitsch         PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T));
1057680d18d5SMatthew G. Knepley         {
1058680d18d5SMatthew G. Knepley           const PetscReal *basis = T->T[0];
1059680d18d5SMatthew G. Knepley           const PetscInt   Nb    = T->Nb;
1060680d18d5SMatthew G. Knepley           const PetscInt   Nc    = T->Nc;
106166a0a883SMatthew G. Knepley           PetscInt         f, fc;
106266a0a883SMatthew G. Knepley 
1063680d18d5SMatthew G. Knepley           for (fc = 0; fc < Nc; ++fc) {
106466a0a883SMatthew G. Knepley             interpolant[p * ctx->dof + coff + fc] = 0.0;
1065ad540459SPierre Jolivet             for (f = 0; f < Nb; ++f) interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc];
1066680d18d5SMatthew G. Knepley           }
106766a0a883SMatthew G. Knepley           coff += Nc;
106866a0a883SMatthew G. Knepley           foff += Nb;
1069680d18d5SMatthew G. Knepley         }
10709566063dSJacob Faibussowitsch         PetscCall(PetscTabulationDestroy(&T));
1071680d18d5SMatthew G. Knepley       }
10729566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
107363a3b9bcSJacob Faibussowitsch       PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof);
107463a3b9bcSJacob Faibussowitsch       PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize);
107566a0a883SMatthew G. Knepley     }
10769566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
10779566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(v, &interpolant));
107866a0a883SMatthew G. Knepley   } else {
107966a0a883SMatthew G. Knepley     DMPolytopeType ct;
108066a0a883SMatthew G. Knepley 
1081680d18d5SMatthew G. Knepley     /* TODO Check each cell individually */
10829566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, ctx->cells[0], &ct));
1083680d18d5SMatthew G. Knepley     switch (ct) {
10849566063dSJacob Faibussowitsch     case DM_POLYTOPE_SEGMENT: PetscCall(DMInterpolate_Segment_Private(ctx, dm, x, v)); break;
10859566063dSJacob Faibussowitsch     case DM_POLYTOPE_TRIANGLE: PetscCall(DMInterpolate_Triangle_Private(ctx, dm, x, v)); break;
10869566063dSJacob Faibussowitsch     case DM_POLYTOPE_QUADRILATERAL: PetscCall(DMInterpolate_Quad_Private(ctx, dm, x, v)); break;
10879566063dSJacob Faibussowitsch     case DM_POLYTOPE_TETRAHEDRON: PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, x, v)); break;
10889566063dSJacob Faibussowitsch     case DM_POLYTOPE_HEXAHEDRON: PetscCall(DMInterpolate_Hex_Private(ctx, dm, x, v)); break;
1089cfe5744fSMatthew G. Knepley     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]);
1090680d18d5SMatthew G. Knepley     }
1091552f7358SJed Brown   }
1092552f7358SJed Brown   PetscFunctionReturn(0);
1093552f7358SJed Brown }
1094552f7358SJed Brown 
10954267b1a3SMatthew G. Knepley /*@C
1096*f6dfbefdSBarry Smith   DMInterpolationDestroy - Destroys a `DMInterpolationInfo` context
10974267b1a3SMatthew G. Knepley 
10984267b1a3SMatthew G. Knepley   Collective on ctx
10994267b1a3SMatthew G. Knepley 
11004267b1a3SMatthew G. Knepley   Input Parameter:
11014267b1a3SMatthew G. Knepley . ctx - the context
11024267b1a3SMatthew G. Knepley 
11034267b1a3SMatthew G. Knepley   Level: beginner
11044267b1a3SMatthew G. Knepley 
1105*f6dfbefdSBarry Smith .seealso: `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
11064267b1a3SMatthew G. Knepley @*/
11079371c9d4SSatish Balay PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) {
1108552f7358SJed Brown   PetscFunctionBegin;
1109064a246eSJacob Faibussowitsch   PetscValidPointer(ctx, 1);
11109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->coords));
11119566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->points));
11129566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->cells));
11139566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
11140298fd71SBarry Smith   *ctx = NULL;
1115552f7358SJed Brown   PetscFunctionReturn(0);
1116552f7358SJed Brown }
1117cc0c4584SMatthew G. Knepley 
1118cc0c4584SMatthew G. Knepley /*@C
1119cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1120cc0c4584SMatthew G. Knepley 
1121*f6dfbefdSBarry Smith   Collective on snes
1122cc0c4584SMatthew G. Knepley 
1123cc0c4584SMatthew G. Knepley   Input Parameters:
1124*f6dfbefdSBarry Smith + snes   - the `SNES` context
1125cc0c4584SMatthew G. Knepley . its    - iteration number
1126cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1127*f6dfbefdSBarry Smith - vf  - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII`
1128cc0c4584SMatthew G. Knepley 
1129*f6dfbefdSBarry Smith   Note:
1130cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1131cc0c4584SMatthew G. Knepley 
1132cc0c4584SMatthew G. Knepley   Level: intermediate
1133cc0c4584SMatthew G. Knepley 
1134*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()`
1135cc0c4584SMatthew G. Knepley @*/
11369371c9d4SSatish Balay PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) {
1137d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1138cc0c4584SMatthew G. Knepley   Vec                res;
1139cc0c4584SMatthew G. Knepley   DM                 dm;
1140cc0c4584SMatthew G. Knepley   PetscSection       s;
1141cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1142cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1143cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1144cc0c4584SMatthew G. Knepley 
1145cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11464d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
11479566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
11489566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
11499566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
11509566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
11519566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
11529566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms));
11539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(res, &r));
1154cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1155cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1156cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1157cc0c4584SMatthew G. Knepley 
11589566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
11599566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
1160cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d]));
1161cc0c4584SMatthew G. Knepley     }
1162cc0c4584SMatthew G. Knepley   }
11639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(res, &r));
11641c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
11659566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
11669566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
116763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm));
1168cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
11699566063dSJacob Faibussowitsch     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
11709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f])));
1171cc0c4584SMatthew G. Knepley   }
11729566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
11739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
11749566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
11759566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lnorms, norms));
1176cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1177cc0c4584SMatthew G. Knepley }
117824cdb843SMatthew G. Knepley 
117924cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
118024cdb843SMatthew G. Knepley 
11819371c9d4SSatish Balay PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS) {
11826528b96dSMatthew G. Knepley   PetscInt depth;
11836528b96dSMatthew G. Knepley 
11846528b96dSMatthew G. Knepley   PetscFunctionBegin;
11859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
11869566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS));
11879566063dSJacob Faibussowitsch   if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS));
11886528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
11896528b96dSMatthew G. Knepley }
11906528b96dSMatthew G. Knepley 
119124cdb843SMatthew G. Knepley /*@
11928db23a53SJed Brown   DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user
119324cdb843SMatthew G. Knepley 
119424cdb843SMatthew G. Knepley   Input Parameters:
119524cdb843SMatthew G. Knepley + dm - The mesh
119624cdb843SMatthew G. Knepley . X  - Local solution
119724cdb843SMatthew G. Knepley - user - The user context
119824cdb843SMatthew G. Knepley 
119924cdb843SMatthew G. Knepley   Output Parameter:
120024cdb843SMatthew G. Knepley . F  - Local output vector
120124cdb843SMatthew G. Knepley 
1202*f6dfbefdSBarry Smith   Note:
1203*f6dfbefdSBarry Smith   The residual is summed into F; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in F is intentional.
12048db23a53SJed Brown 
120524cdb843SMatthew G. Knepley   Level: developer
120624cdb843SMatthew G. Knepley 
1207*f6dfbefdSBarry Smith .seealso: `DM`, `DMPlexComputeJacobianAction()`
120824cdb843SMatthew G. Knepley @*/
12099371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) {
12106da023fcSToby Isaac   DM       plex;
1211083401c6SMatthew G. Knepley   IS       allcellIS;
12126528b96dSMatthew G. Knepley   PetscInt Nds, s;
121324cdb843SMatthew G. Knepley 
121424cdb843SMatthew G. Knepley   PetscFunctionBegin;
12159566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12169566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12179566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
12186528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
12196528b96dSMatthew G. Knepley     PetscDS      ds;
12206528b96dSMatthew G. Knepley     IS           cellIS;
122106ad1575SMatthew G. Knepley     PetscFormKey key;
12226528b96dSMatthew G. Knepley 
12239566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
12246528b96dSMatthew G. Knepley     key.value = 0;
12256528b96dSMatthew G. Knepley     key.field = 0;
122606ad1575SMatthew G. Knepley     key.part  = 0;
12276528b96dSMatthew G. Knepley     if (!key.label) {
12289566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
12296528b96dSMatthew G. Knepley       cellIS = allcellIS;
12306528b96dSMatthew G. Knepley     } else {
12316528b96dSMatthew G. Knepley       IS pointIS;
12326528b96dSMatthew G. Knepley 
12336528b96dSMatthew G. Knepley       key.value = 1;
12349566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
12359566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12369566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
12376528b96dSMatthew G. Knepley     }
12389566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
12399566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
12406528b96dSMatthew G. Knepley   }
12419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
12429566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
12436528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
12446528b96dSMatthew G. Knepley }
12456528b96dSMatthew G. Knepley 
12469371c9d4SSatish Balay PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user) {
12476528b96dSMatthew G. Knepley   DM       plex;
12486528b96dSMatthew G. Knepley   IS       allcellIS;
12496528b96dSMatthew G. Knepley   PetscInt Nds, s;
12506528b96dSMatthew G. Knepley 
12516528b96dSMatthew G. Knepley   PetscFunctionBegin;
12529566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12539566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12549566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1255083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1256083401c6SMatthew G. Knepley     PetscDS ds;
1257083401c6SMatthew G. Knepley     DMLabel label;
1258083401c6SMatthew G. Knepley     IS      cellIS;
1259083401c6SMatthew G. Knepley 
12609566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
12616528b96dSMatthew G. Knepley     {
126206ad1575SMatthew G. Knepley       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
12636528b96dSMatthew G. Knepley       PetscWeakForm     wf;
12646528b96dSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
126506ad1575SMatthew G. Knepley       PetscFormKey     *reskeys;
12666528b96dSMatthew G. Knepley 
12676528b96dSMatthew G. Knepley       /* Get unique residual keys */
12686528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
12696528b96dSMatthew G. Knepley         PetscInt Nkm;
12709566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm));
12716528b96dSMatthew G. Knepley         Nk += Nkm;
12726528b96dSMatthew G. Knepley       }
12739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &reskeys));
127448a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys));
127563a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
12769566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, reskeys));
12776528b96dSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
12786528b96dSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
12796528b96dSMatthew G. Knepley           ++k;
12806528b96dSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
12816528b96dSMatthew G. Knepley         }
12826528b96dSMatthew G. Knepley       }
12836528b96dSMatthew G. Knepley       Nk = k;
12846528b96dSMatthew G. Knepley 
12859566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
12866528b96dSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
12876528b96dSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
12886528b96dSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
12896528b96dSMatthew G. Knepley 
1290083401c6SMatthew G. Knepley         if (!label) {
12919566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
1292083401c6SMatthew G. Knepley           cellIS = allcellIS;
1293083401c6SMatthew G. Knepley         } else {
1294083401c6SMatthew G. Knepley           IS pointIS;
1295083401c6SMatthew G. Knepley 
12969566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
12979566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12989566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
12994a3e9fdbSToby Isaac         }
13009566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
13019566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
1302083401c6SMatthew G. Knepley       }
13039566063dSJacob Faibussowitsch       PetscCall(PetscFree(reskeys));
13046528b96dSMatthew G. Knepley     }
13056528b96dSMatthew G. Knepley   }
13069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
13079566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
130824cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
130924cdb843SMatthew G. Knepley }
131024cdb843SMatthew G. Knepley 
1311bdd6f66aSToby Isaac /*@
1312bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1313bdd6f66aSToby Isaac 
1314bdd6f66aSToby Isaac   Input Parameters:
1315bdd6f66aSToby Isaac + dm - The mesh
1316bdd6f66aSToby Isaac - user - The user context
1317bdd6f66aSToby Isaac 
1318bdd6f66aSToby Isaac   Output Parameter:
1319bdd6f66aSToby Isaac . X  - Local solution
1320bdd6f66aSToby Isaac 
1321bdd6f66aSToby Isaac   Level: developer
1322bdd6f66aSToby Isaac 
1323*f6dfbefdSBarry Smith .seealso: `DMPLEX`, `DMPlexComputeJacobianAction()`
1324bdd6f66aSToby Isaac @*/
13259371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) {
1326bdd6f66aSToby Isaac   DM plex;
1327bdd6f66aSToby Isaac 
1328bdd6f66aSToby Isaac   PetscFunctionBegin;
13299566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
13309566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL));
13319566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1332bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1333bdd6f66aSToby Isaac }
1334bdd6f66aSToby Isaac 
13357a73cf09SMatthew G. Knepley /*@
13368e3a2eefSMatthew G. Knepley   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
13377a73cf09SMatthew G. Knepley 
13387a73cf09SMatthew G. Knepley   Input Parameters:
1339*f6dfbefdSBarry Smith + dm   - The `DM`
13407a73cf09SMatthew G. Knepley . X    - Local solution vector
13417a73cf09SMatthew G. Knepley . Y    - Local input vector
13427a73cf09SMatthew G. Knepley - user - The user context
13437a73cf09SMatthew G. Knepley 
13447a73cf09SMatthew G. Knepley   Output Parameter:
13458e3a2eefSMatthew G. Knepley . F    - lcoal output vector
13467a73cf09SMatthew G. Knepley 
13477a73cf09SMatthew G. Knepley   Level: developer
13487a73cf09SMatthew G. Knepley 
13498e3a2eefSMatthew G. Knepley   Notes:
1350*f6dfbefdSBarry Smith   Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly.
13518e3a2eefSMatthew G. Knepley 
1352*f6dfbefdSBarry Smith .seealso: `DM`, ``DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()`
13537a73cf09SMatthew G. Knepley @*/
13549371c9d4SSatish Balay PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user) {
13558e3a2eefSMatthew G. Knepley   DM       plex;
13568e3a2eefSMatthew G. Knepley   IS       allcellIS;
13578e3a2eefSMatthew G. Knepley   PetscInt Nds, s;
1358a925c78cSMatthew G. Knepley 
1359a925c78cSMatthew G. Knepley   PetscFunctionBegin;
13609566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
13619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
13629566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
13638e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
13648e3a2eefSMatthew G. Knepley     PetscDS ds;
13658e3a2eefSMatthew G. Knepley     DMLabel label;
13668e3a2eefSMatthew G. Knepley     IS      cellIS;
13677a73cf09SMatthew G. Knepley 
13689566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
13698e3a2eefSMatthew G. Knepley     {
137006ad1575SMatthew G. Knepley       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
13718e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
13728e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
137306ad1575SMatthew G. Knepley       PetscFormKey     *jackeys;
13748e3a2eefSMatthew G. Knepley 
13758e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
13768e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
13778e3a2eefSMatthew G. Knepley         PetscInt Nkm;
13789566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm));
13798e3a2eefSMatthew G. Knepley         Nk += Nkm;
13808e3a2eefSMatthew G. Knepley       }
13819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &jackeys));
138248a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys));
138363a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
13849566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, jackeys));
13858e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
13868e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
13878e3a2eefSMatthew G. Knepley           ++k;
13888e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
13898e3a2eefSMatthew G. Knepley         }
13908e3a2eefSMatthew G. Knepley       }
13918e3a2eefSMatthew G. Knepley       Nk = k;
13928e3a2eefSMatthew G. Knepley 
13939566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
13948e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
13958e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
13968e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
13978e3a2eefSMatthew G. Knepley 
13988e3a2eefSMatthew G. Knepley         if (!label) {
13999566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
14008e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
14017a73cf09SMatthew G. Knepley         } else {
14028e3a2eefSMatthew G. Knepley           IS pointIS;
1403a925c78cSMatthew G. Knepley 
14049566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
14059566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14069566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
1407a925c78cSMatthew G. Knepley         }
14089566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user));
14099566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
14108e3a2eefSMatthew G. Knepley       }
14119566063dSJacob Faibussowitsch       PetscCall(PetscFree(jackeys));
14128e3a2eefSMatthew G. Knepley     }
14138e3a2eefSMatthew G. Knepley   }
14149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
14159566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1416a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1417a925c78cSMatthew G. Knepley }
1418a925c78cSMatthew G. Knepley 
141924cdb843SMatthew G. Knepley /*@
142024cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
142124cdb843SMatthew G. Knepley 
142224cdb843SMatthew G. Knepley   Input Parameters:
142324cdb843SMatthew G. Knepley + dm - The mesh
142424cdb843SMatthew G. Knepley . X  - Local input vector
142524cdb843SMatthew G. Knepley - user - The user context
142624cdb843SMatthew G. Knepley 
142724cdb843SMatthew G. Knepley   Output Parameter:
142824cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
142924cdb843SMatthew G. Knepley 
143024cdb843SMatthew G. Knepley   Note:
143124cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
143224cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
143324cdb843SMatthew G. Knepley 
143424cdb843SMatthew G. Knepley   Level: developer
143524cdb843SMatthew G. Knepley 
1436*f6dfbefdSBarry Smith .seealso: `DMPLEX`, `Mat`
143724cdb843SMatthew G. Knepley @*/
14389371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user) {
14396da023fcSToby Isaac   DM        plex;
1440083401c6SMatthew G. Knepley   IS        allcellIS;
1441f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
14426528b96dSMatthew G. Knepley   PetscInt  Nds, s;
144324cdb843SMatthew G. Knepley 
144424cdb843SMatthew G. Knepley   PetscFunctionBegin;
14459566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
14469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
14479566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1448083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1449083401c6SMatthew G. Knepley     PetscDS      ds;
1450083401c6SMatthew G. Knepley     IS           cellIS;
145106ad1575SMatthew G. Knepley     PetscFormKey key;
1452083401c6SMatthew G. Knepley 
14539566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
14546528b96dSMatthew G. Knepley     key.value = 0;
14556528b96dSMatthew G. Knepley     key.field = 0;
145606ad1575SMatthew G. Knepley     key.part  = 0;
14576528b96dSMatthew G. Knepley     if (!key.label) {
14589566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1459083401c6SMatthew G. Knepley       cellIS = allcellIS;
1460083401c6SMatthew G. Knepley     } else {
1461083401c6SMatthew G. Knepley       IS pointIS;
1462083401c6SMatthew G. Knepley 
14636528b96dSMatthew G. Knepley       key.value = 1;
14649566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
14659566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14669566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1467083401c6SMatthew G. Knepley     }
1468083401c6SMatthew G. Knepley     if (!s) {
14699566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
14709566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
14719566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
14729566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
1473083401c6SMatthew G. Knepley     }
14749566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user));
14759566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1476083401c6SMatthew G. Knepley   }
14779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
14789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
147924cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
148024cdb843SMatthew G. Knepley }
14811878804aSMatthew G. Knepley 
14829371c9d4SSatish Balay struct _DMSNESJacobianMFCtx {
14838e3a2eefSMatthew G. Knepley   DM    dm;
14848e3a2eefSMatthew G. Knepley   Vec   X;
14858e3a2eefSMatthew G. Knepley   void *ctx;
14868e3a2eefSMatthew G. Knepley };
14878e3a2eefSMatthew G. Knepley 
14889371c9d4SSatish Balay static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A) {
14898e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
14908e3a2eefSMatthew G. Knepley 
14918e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
14929566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
14939566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, NULL));
14949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
14959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->X));
14969566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
14978e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
14988e3a2eefSMatthew G. Knepley }
14998e3a2eefSMatthew G. Knepley 
15009371c9d4SSatish Balay static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z) {
15018e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15028e3a2eefSMatthew G. Knepley 
15038e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15049566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
15059566063dSJacob Faibussowitsch   PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx));
15068e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15078e3a2eefSMatthew G. Knepley }
15088e3a2eefSMatthew G. Knepley 
15098e3a2eefSMatthew G. Knepley /*@
1510*f6dfbefdSBarry Smith   DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free
15118e3a2eefSMatthew G. Knepley 
15128e3a2eefSMatthew G. Knepley   Collective on dm
15138e3a2eefSMatthew G. Knepley 
15148e3a2eefSMatthew G. Knepley   Input Parameters:
1515*f6dfbefdSBarry Smith + dm   - The `DM`
15168e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
15178e3a2eefSMatthew G. Knepley - user - A user context, or NULL
15188e3a2eefSMatthew G. Knepley 
15198e3a2eefSMatthew G. Knepley   Output Parameter:
1520*f6dfbefdSBarry Smith . J    - The `Mat`
15218e3a2eefSMatthew G. Knepley 
15228e3a2eefSMatthew G. Knepley   Level: advanced
15238e3a2eefSMatthew G. Knepley 
1524*f6dfbefdSBarry Smith   Note:
1525*f6dfbefdSBarry Smith   Vec X is kept in `Mat` J, so updating X then updates the evaluation point.
15268e3a2eefSMatthew G. Knepley 
1527*f6dfbefdSBarry Smith .seealso: `DM`, `DMSNESComputeJacobianAction()`
15288e3a2eefSMatthew G. Knepley @*/
15299371c9d4SSatish Balay PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J) {
15308e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15318e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
15328e3a2eefSMatthew G. Knepley 
15338e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15349566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
15359566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
15369566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
15379566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
15389566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, n, n, N, N));
15399566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
15409566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)X));
15419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &ctx));
15428e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
15438e3a2eefSMatthew G. Knepley   ctx->X   = X;
15448e3a2eefSMatthew G. Knepley   ctx->ctx = user;
15459566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, ctx));
15469566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private));
15479566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private));
15488e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15498e3a2eefSMatthew G. Knepley }
15508e3a2eefSMatthew G. Knepley 
155138cfc46eSPierre Jolivet /*
155238cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
155338cfc46eSPierre Jolivet 
155438cfc46eSPierre Jolivet    Input Parameters:
155538cfc46eSPierre Jolivet +     X - SNES linearization point
155638cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
155738cfc46eSPierre Jolivet 
155838cfc46eSPierre Jolivet    Output Parameter:
155938cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
156038cfc46eSPierre Jolivet 
156138cfc46eSPierre Jolivet    Level: intermediate
156238cfc46eSPierre Jolivet 
1563db781477SPatrick Sanan .seealso: `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()`
156438cfc46eSPierre Jolivet */
15659371c9d4SSatish Balay static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) {
156638cfc46eSPierre Jolivet   SNES   snes;
156738cfc46eSPierre Jolivet   Mat    pJ;
156838cfc46eSPierre Jolivet   DM     ovldm, origdm;
156938cfc46eSPierre Jolivet   DMSNES sdm;
157038cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM, Vec, void *);
157138cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *);
157238cfc46eSPierre Jolivet   void *bctx, *jctx;
157338cfc46eSPierre Jolivet 
157438cfc46eSPierre Jolivet   PetscFunctionBegin;
15759566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ));
157628b400f6SJacob Faibussowitsch   PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat");
15779566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm));
157828b400f6SJacob Faibussowitsch   PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM");
15799566063dSJacob Faibussowitsch   PetscCall(MatGetDM(pJ, &ovldm));
15809566063dSJacob Faibussowitsch   PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx));
15819566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx));
15829566063dSJacob Faibussowitsch   PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx));
15839566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx));
15849566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes));
158538cfc46eSPierre Jolivet   if (!snes) {
15869566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes));
15879566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, ovldm));
15889566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes));
15899566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)snes));
159038cfc46eSPierre Jolivet   }
15919566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(ovldm, &sdm));
15929566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(X));
1593800f99ffSJeremy L Thompson   {
1594800f99ffSJeremy L Thompson     void *ctx;
1595800f99ffSJeremy L Thompson     PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *);
1596800f99ffSJeremy L Thompson     PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx));
1597800f99ffSJeremy L Thompson     PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx));
1598800f99ffSJeremy L Thompson   }
15999566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(X));
160038cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
160138cfc46eSPierre Jolivet   {
160238cfc46eSPierre Jolivet     Mat locpJ;
160338cfc46eSPierre Jolivet 
16049566063dSJacob Faibussowitsch     PetscCall(MatISGetLocalMat(pJ, &locpJ));
16059566063dSJacob Faibussowitsch     PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN));
160638cfc46eSPierre Jolivet   }
160738cfc46eSPierre Jolivet   PetscFunctionReturn(0);
160838cfc46eSPierre Jolivet }
160938cfc46eSPierre Jolivet 
1610a925c78cSMatthew G. Knepley /*@
1611*f6dfbefdSBarry Smith   DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, residual, and Jacobian.
16129f520fc2SToby Isaac 
16139f520fc2SToby Isaac   Input Parameters:
1614*f6dfbefdSBarry Smith + dm - The `DM` object
1615*f6dfbefdSBarry Smith . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see `PetscDSAddBoundary()`)
1616*f6dfbefdSBarry Smith . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see `PetscDSSetResidual()`)
1617*f6dfbefdSBarry Smith - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see `PetscDSSetJacobian()`)
16181a244344SSatish Balay 
16191a244344SSatish Balay   Level: developer
1620*f6dfbefdSBarry Smith 
1621*f6dfbefdSBarry Smith .seealso: `DMPLEX`, `SNES`
16229f520fc2SToby Isaac @*/
16239371c9d4SSatish Balay PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) {
16249f520fc2SToby Isaac   PetscFunctionBegin;
16259566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, boundaryctx));
16269566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, residualctx));
16279566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, jacobianctx));
16289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex));
16299f520fc2SToby Isaac   PetscFunctionReturn(0);
16309f520fc2SToby Isaac }
16319f520fc2SToby Isaac 
1632f017bcb6SMatthew G. Knepley /*@C
1633f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1634f017bcb6SMatthew G. Knepley 
1635f017bcb6SMatthew G. Knepley   Input Parameters:
1636*f6dfbefdSBarry Smith + snes - the `SNES` object
1637*f6dfbefdSBarry Smith . dm   - the `DM`
1638f2cacb80SMatthew G. Knepley . t    - the time
1639*f6dfbefdSBarry Smith . u    - a `DM` vector
1640f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1641f017bcb6SMatthew G. Knepley 
1642f3c94560SMatthew G. Knepley   Output Parameters:
1643f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1644f3c94560SMatthew G. Knepley 
1645*f6dfbefdSBarry Smith   Note:
1646*f6dfbefdSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
16477f96f943SMatthew G. Knepley 
1648f017bcb6SMatthew G. Knepley   Level: developer
1649f017bcb6SMatthew G. Knepley 
1650*f6dfbefdSBarry Smith .seealso: `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`, `PetscDSSetExactSolution()`
1651f017bcb6SMatthew G. Knepley @*/
16529371c9d4SSatish Balay PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[]) {
1653f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1654f017bcb6SMatthew G. Knepley   void     **ectxs;
1655f3c94560SMatthew G. Knepley   PetscReal *err;
16567f96f943SMatthew G. Knepley   MPI_Comm   comm;
16577f96f943SMatthew G. Knepley   PetscInt   Nf, f;
16581878804aSMatthew G. Knepley 
16591878804aSMatthew G. Knepley   PetscFunctionBegin;
1660f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1661f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1662064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1663534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
16647f96f943SMatthew G. Knepley 
16659566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
16669566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
16677f96f943SMatthew G. Knepley 
16689566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
16699566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
16709566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
16717f96f943SMatthew G. Knepley   {
16727f96f943SMatthew G. Knepley     PetscInt Nds, s;
16737f96f943SMatthew G. Knepley 
16749566063dSJacob Faibussowitsch     PetscCall(DMGetNumDS(dm, &Nds));
1675083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
16767f96f943SMatthew G. Knepley       PetscDS         ds;
1677083401c6SMatthew G. Knepley       DMLabel         label;
1678083401c6SMatthew G. Knepley       IS              fieldIS;
16797f96f943SMatthew G. Knepley       const PetscInt *fields;
16807f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1681083401c6SMatthew G. Knepley 
16829566063dSJacob Faibussowitsch       PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds));
16839566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(ds, &dsNf));
16849566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
1685083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1686083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
16879566063dSJacob Faibussowitsch         PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
1688083401c6SMatthew G. Knepley       }
16899566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
1690083401c6SMatthew G. Knepley     }
1691083401c6SMatthew G. Knepley   }
1692f017bcb6SMatthew G. Knepley   if (Nf > 1) {
16939566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err));
1694f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1695ad540459SPierre 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);
1696b878532fSMatthew G. Knepley     } else if (error) {
1697b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
16981878804aSMatthew G. Knepley     } else {
16999566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: ["));
1700f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
17019566063dSJacob Faibussowitsch         if (f) PetscCall(PetscPrintf(comm, ", "));
17029566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(comm, "%g", (double)err[f]));
17031878804aSMatthew G. Knepley       }
17049566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "]\n"));
1705f017bcb6SMatthew G. Knepley     }
1706f017bcb6SMatthew G. Knepley   } else {
17079566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]));
1708f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
170908401ef6SPierre Jolivet       PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol);
1710b878532fSMatthew G. Knepley     } else if (error) {
1711b878532fSMatthew G. Knepley       error[0] = err[0];
1712f017bcb6SMatthew G. Knepley     } else {
17139566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]));
1714f017bcb6SMatthew G. Knepley     }
1715f017bcb6SMatthew G. Knepley   }
17169566063dSJacob Faibussowitsch   PetscCall(PetscFree3(exacts, ectxs, err));
1717f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1718f017bcb6SMatthew G. Knepley }
1719f017bcb6SMatthew G. Knepley 
1720f017bcb6SMatthew G. Knepley /*@C
1721f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1722f017bcb6SMatthew G. Knepley 
1723f017bcb6SMatthew G. Knepley   Input Parameters:
1724*f6dfbefdSBarry Smith + snes - the `SNES` object
1725*f6dfbefdSBarry Smith . dm   - the `DM`
1726*f6dfbefdSBarry Smith . u    - a `DM` vector
1727f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1728f017bcb6SMatthew G. Knepley 
1729*f6dfbefdSBarry Smith   Output Parameter:
1730f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1731f3c94560SMatthew G. Knepley 
1732f017bcb6SMatthew G. Knepley   Level: developer
1733f017bcb6SMatthew G. Knepley 
1734db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
1735f017bcb6SMatthew G. Knepley @*/
17369371c9d4SSatish Balay PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) {
1737f017bcb6SMatthew G. Knepley   MPI_Comm  comm;
1738f017bcb6SMatthew G. Knepley   Vec       r;
1739f017bcb6SMatthew G. Knepley   PetscReal res;
1740f017bcb6SMatthew G. Knepley 
1741f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1742f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1743f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1744f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1745534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
17469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
17479566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
17489566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
17499566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, u, r));
17509566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
1751f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
175208401ef6SPierre Jolivet     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
1753b878532fSMatthew G. Knepley   } else if (residual) {
1754b878532fSMatthew G. Knepley     *residual = res;
1755f017bcb6SMatthew G. Knepley   } else {
17569566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
17579566063dSJacob Faibussowitsch     PetscCall(VecChop(r, 1.0e-10));
17589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
17599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
17609566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
1761f017bcb6SMatthew G. Knepley   }
17629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1763f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1764f017bcb6SMatthew G. Knepley }
1765f017bcb6SMatthew G. Knepley 
1766f017bcb6SMatthew G. Knepley /*@C
1767f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1768f017bcb6SMatthew G. Knepley 
1769f017bcb6SMatthew G. Knepley   Input Parameters:
1770*f6dfbefdSBarry Smith + snes - the `SNES` object
1771*f6dfbefdSBarry Smith . dm   - the `DM`
1772*f6dfbefdSBarry Smith . u    - a `DM` vector
1773f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1774f017bcb6SMatthew G. Knepley 
1775f3c94560SMatthew G. Knepley   Output Parameters:
1776f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1777f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1778f3c94560SMatthew G. Knepley 
1779f017bcb6SMatthew G. Knepley   Level: developer
1780f017bcb6SMatthew G. Knepley 
1781db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
1782f017bcb6SMatthew G. Knepley @*/
17839371c9d4SSatish Balay PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) {
1784f017bcb6SMatthew G. Knepley   MPI_Comm     comm;
1785f017bcb6SMatthew G. Knepley   PetscDS      ds;
1786f017bcb6SMatthew G. Knepley   Mat          J, M;
1787f017bcb6SMatthew G. Knepley   MatNullSpace nullspace;
1788f3c94560SMatthew G. Knepley   PetscReal    slope, intercept;
17897f96f943SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
1790f017bcb6SMatthew G. Knepley 
1791f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1792f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1793f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1794f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1795534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1796064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 6);
17979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
17989566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
1799f017bcb6SMatthew G. Knepley   /* Create and view matrices */
18009566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
18019566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
18029566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
18039566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
1804282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
18059566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
18069566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, M));
18079566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
18089566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
18099566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
18109566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
1811282e7bb4SMatthew G. Knepley   } else {
18129566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, J));
1813282e7bb4SMatthew G. Knepley   }
18149566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
18159566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
18169566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
1817f017bcb6SMatthew G. Knepley   /* Check nullspace */
18189566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
1819f017bcb6SMatthew G. Knepley   if (nullspace) {
18201878804aSMatthew G. Knepley     PetscBool isNull;
18219566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
182228b400f6SJacob Faibussowitsch     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
18231878804aSMatthew G. Knepley   }
1824f017bcb6SMatthew G. Knepley   /* Taylor test */
1825f017bcb6SMatthew G. Knepley   {
1826f017bcb6SMatthew G. Knepley     PetscRandom rand;
1827f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1828f3c94560SMatthew G. Knepley     PetscReal   h;
1829f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1830f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1831f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1832f017bcb6SMatthew G. Knepley 
1833f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
18349566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
18359566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
18369566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
18379566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
18389566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
18399566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
1840f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
18419566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
18429566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, u, r));
1843f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
18449371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
18459371c9d4SSatish Balay       ;
18469566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
18479566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
18489566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
1849f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
18509566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
1851f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
18529566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, uhat, rhat));
18539566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
18549566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
1855f017bcb6SMatthew G. Knepley 
1856f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1857f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1858f017bcb6SMatthew G. Knepley     }
18599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
18609566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
18619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
18629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
18639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
1864f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1865f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1866f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
1867f017bcb6SMatthew G. Knepley     }
1868f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
18699566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
18709566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
1871f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1872f017bcb6SMatthew G. Knepley     if (tol >= 0) {
18730b121fc5SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
1874b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1875b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1876b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1877f017bcb6SMatthew G. Knepley     } else {
18789566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
18799566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
1880f017bcb6SMatthew G. Knepley     }
1881f017bcb6SMatthew G. Knepley   }
18829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
18831878804aSMatthew G. Knepley   PetscFunctionReturn(0);
18841878804aSMatthew G. Knepley }
18851878804aSMatthew G. Knepley 
18869371c9d4SSatish Balay PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u) {
1887f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
18889566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL));
18899566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL));
18909566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL));
1891f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1892f017bcb6SMatthew G. Knepley }
1893f017bcb6SMatthew G. Knepley 
1894bee9a294SMatthew G. Knepley /*@C
1895bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1896bee9a294SMatthew G. Knepley 
1897bee9a294SMatthew G. Knepley   Input Parameters:
1898*f6dfbefdSBarry Smith + snes - the `SNES` object
1899*f6dfbefdSBarry Smith - u    - representative `SNES` vector
19007f96f943SMatthew G. Knepley 
1901*f6dfbefdSBarry Smith   Note:
1902*f6dfbefdSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
1903bee9a294SMatthew G. Knepley 
1904bee9a294SMatthew G. Knepley   Level: developer
1905bee9a294SMatthew G. Knepley @*/
19069371c9d4SSatish Balay PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u) {
19071878804aSMatthew G. Knepley   DM        dm;
19081878804aSMatthew G. Knepley   Vec       sol;
19091878804aSMatthew G. Knepley   PetscBool check;
19101878804aSMatthew G. Knepley 
19111878804aSMatthew G. Knepley   PetscFunctionBegin;
19129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check));
19131878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
19149566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
19159566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
19169566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, sol));
19179566063dSJacob Faibussowitsch   PetscCall(DMSNESCheck_Internal(snes, dm, sol));
19189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
1919552f7358SJed Brown   PetscFunctionReturn(0);
1920552f7358SJed Brown }
1921