xref: /petsc/src/snes/utils/dmplexsnes.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
14649ef022SMatthew Knepley   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
66649ef022SMatthew Knepley    SNESConvergedCorrectPressure - Convergence test that adds 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. The convergence test itself just mimics SNESConvergedDefault().
67649ef022SMatthew Knepley 
68649ef022SMatthew Knepley    Logically Collective on SNES
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley    Input Parameters:
71649ef022SMatthew Knepley +  snes - the SNES context
72649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
73649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
74649ef022SMatthew Knepley .  snorm - 2-norm of current step
75649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
76649ef022SMatthew Knepley -  ctx   - Optional user context
77649ef022SMatthew Knepley 
78649ef022SMatthew Knepley    Output Parameter:
79649ef022SMatthew Knepley .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
80649ef022SMatthew Knepley 
81649ef022SMatthew Knepley    Notes:
82f362779dSJed Brown    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 must be created with discretizations of those fields. We currently assume that the pressure field has index 1. The pressure field must have a nullspace, likely created using the DMSetNullSpaceConstructor() interface. Last we must be able to integrate the pressure over the domain, so the DM attached to the SNES must be a Plex at this time.
83f362779dSJed Brown 
84f362779dSJed Brown    Options Database Keys:
85f362779dSJed Brown .  -snes_convergence_test correct_pressure - see SNESSetFromOptions()
86649ef022SMatthew Knepley 
87649ef022SMatthew Knepley    Level: advanced
88649ef022SMatthew Knepley 
89db781477SPatrick Sanan .seealso: `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()`
90649ef022SMatthew Knepley @*/
919371c9d4SSatish Balay PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx) {
92649ef022SMatthew Knepley   PetscBool monitorIntegral = PETSC_FALSE;
93649ef022SMatthew Knepley 
94649ef022SMatthew Knepley   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx));
96649ef022SMatthew Knepley   if (monitorIntegral) {
97649ef022SMatthew Knepley     Mat          J;
98649ef022SMatthew Knepley     Vec          u;
99649ef022SMatthew Knepley     MatNullSpace nullspace;
100649ef022SMatthew Knepley     const Vec   *nullvecs;
101649ef022SMatthew Knepley     PetscScalar  pintd;
102649ef022SMatthew Knepley 
1039566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1049566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1059566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1069566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs));
1079566063dSJacob Faibussowitsch     PetscCall(VecDot(nullvecs[0], u, &pintd));
1089566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double)PetscRealPart(pintd)));
109649ef022SMatthew Knepley   }
110649ef022SMatthew Knepley   if (*reason > 0) {
111649ef022SMatthew Knepley     Mat          J;
112649ef022SMatthew Knepley     Vec          u;
113649ef022SMatthew Knepley     MatNullSpace nullspace;
114649ef022SMatthew Knepley     PetscInt     pfield = 1;
115649ef022SMatthew Knepley 
1169566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1179566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1189566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1199566063dSJacob Faibussowitsch     PetscCall(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx));
120649ef022SMatthew Knepley   }
121649ef022SMatthew Knepley   PetscFunctionReturn(0);
122649ef022SMatthew Knepley }
123649ef022SMatthew Knepley 
12424cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
12524cdb843SMatthew G. Knepley 
1269371c9d4SSatish Balay static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) {
1276da023fcSToby Isaac   PetscBool isPlex;
1286da023fcSToby Isaac 
1296da023fcSToby Isaac   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
1316da023fcSToby Isaac   if (isPlex) {
1326da023fcSToby Isaac     *plex = dm;
1339566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
134f7148743SMatthew G. Knepley   } else {
1359566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
136f7148743SMatthew G. Knepley     if (!*plex) {
1379566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
1389566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
1396da023fcSToby Isaac       if (copy) {
1409566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
1419566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
1426da023fcSToby Isaac       }
143f7148743SMatthew G. Knepley     } else {
1449566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*plex));
145f7148743SMatthew G. Knepley     }
1466da023fcSToby Isaac   }
1476da023fcSToby Isaac   PetscFunctionReturn(0);
1486da023fcSToby Isaac }
1496da023fcSToby Isaac 
1504267b1a3SMatthew G. Knepley /*@C
1514267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
1524267b1a3SMatthew G. Knepley 
153d083f849SBarry Smith   Collective
1544267b1a3SMatthew G. Knepley 
1554267b1a3SMatthew G. Knepley   Input Parameter:
1564267b1a3SMatthew G. Knepley . comm - the communicator
1574267b1a3SMatthew G. Knepley 
1584267b1a3SMatthew G. Knepley   Output Parameter:
1594267b1a3SMatthew G. Knepley . ctx - the context
1604267b1a3SMatthew G. Knepley 
1614267b1a3SMatthew G. Knepley   Level: beginner
1624267b1a3SMatthew G. Knepley 
163db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()`
1644267b1a3SMatthew G. Knepley @*/
1659371c9d4SSatish Balay PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) {
166552f7358SJed Brown   PetscFunctionBegin;
167552f7358SJed Brown   PetscValidPointer(ctx, 2);
1689566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
1691aa26658SKarl Rupp 
170552f7358SJed Brown   (*ctx)->comm   = comm;
171552f7358SJed Brown   (*ctx)->dim    = -1;
172552f7358SJed Brown   (*ctx)->nInput = 0;
1730298fd71SBarry Smith   (*ctx)->points = NULL;
1740298fd71SBarry Smith   (*ctx)->cells  = NULL;
175552f7358SJed Brown   (*ctx)->n      = -1;
1760298fd71SBarry Smith   (*ctx)->coords = NULL;
177552f7358SJed Brown   PetscFunctionReturn(0);
178552f7358SJed Brown }
179552f7358SJed Brown 
1804267b1a3SMatthew G. Knepley /*@C
1814267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
1824267b1a3SMatthew G. Knepley 
1834267b1a3SMatthew G. Knepley   Not collective
1844267b1a3SMatthew G. Knepley 
1854267b1a3SMatthew G. Knepley   Input Parameters:
1864267b1a3SMatthew G. Knepley + ctx - the context
1874267b1a3SMatthew G. Knepley - dim - the spatial dimension
1884267b1a3SMatthew G. Knepley 
1894267b1a3SMatthew G. Knepley   Level: intermediate
1904267b1a3SMatthew G. Knepley 
191db781477SPatrick Sanan .seealso: `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
1924267b1a3SMatthew G. Knepley @*/
1939371c9d4SSatish Balay PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) {
194552f7358SJed Brown   PetscFunctionBegin;
19563a3b9bcSJacob Faibussowitsch   PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim);
196552f7358SJed Brown   ctx->dim = dim;
197552f7358SJed Brown   PetscFunctionReturn(0);
198552f7358SJed Brown }
199552f7358SJed Brown 
2004267b1a3SMatthew G. Knepley /*@C
2014267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2024267b1a3SMatthew G. Knepley 
2034267b1a3SMatthew G. Knepley   Not collective
2044267b1a3SMatthew G. Knepley 
2054267b1a3SMatthew G. Knepley   Input Parameter:
2064267b1a3SMatthew G. Knepley . ctx - the context
2074267b1a3SMatthew G. Knepley 
2084267b1a3SMatthew G. Knepley   Output Parameter:
2094267b1a3SMatthew G. Knepley . dim - the spatial dimension
2104267b1a3SMatthew G. Knepley 
2114267b1a3SMatthew G. Knepley   Level: intermediate
2124267b1a3SMatthew G. Knepley 
213db781477SPatrick Sanan .seealso: `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2144267b1a3SMatthew G. Knepley @*/
2159371c9d4SSatish Balay PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) {
216552f7358SJed Brown   PetscFunctionBegin;
217552f7358SJed Brown   PetscValidIntPointer(dim, 2);
218552f7358SJed Brown   *dim = ctx->dim;
219552f7358SJed Brown   PetscFunctionReturn(0);
220552f7358SJed Brown }
221552f7358SJed Brown 
2224267b1a3SMatthew G. Knepley /*@C
2234267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2244267b1a3SMatthew G. Knepley 
2254267b1a3SMatthew G. Knepley   Not collective
2264267b1a3SMatthew G. Knepley 
2274267b1a3SMatthew G. Knepley   Input Parameters:
2284267b1a3SMatthew G. Knepley + ctx - the context
2294267b1a3SMatthew G. Knepley - dof - the number of fields
2304267b1a3SMatthew G. Knepley 
2314267b1a3SMatthew G. Knepley   Level: intermediate
2324267b1a3SMatthew G. Knepley 
233db781477SPatrick Sanan .seealso: `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2344267b1a3SMatthew G. Knepley @*/
2359371c9d4SSatish Balay PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) {
236552f7358SJed Brown   PetscFunctionBegin;
23763a3b9bcSJacob Faibussowitsch   PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof);
238552f7358SJed Brown   ctx->dof = dof;
239552f7358SJed Brown   PetscFunctionReturn(0);
240552f7358SJed Brown }
241552f7358SJed Brown 
2424267b1a3SMatthew G. Knepley /*@C
2434267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2444267b1a3SMatthew G. Knepley 
2454267b1a3SMatthew G. Knepley   Not collective
2464267b1a3SMatthew G. Knepley 
2474267b1a3SMatthew G. Knepley   Input Parameter:
2484267b1a3SMatthew G. Knepley . ctx - the context
2494267b1a3SMatthew G. Knepley 
2504267b1a3SMatthew G. Knepley   Output Parameter:
2514267b1a3SMatthew G. Knepley . dof - the number of fields
2524267b1a3SMatthew G. Knepley 
2534267b1a3SMatthew G. Knepley   Level: intermediate
2544267b1a3SMatthew G. Knepley 
255db781477SPatrick Sanan .seealso: `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2564267b1a3SMatthew G. Knepley @*/
2579371c9d4SSatish Balay PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) {
258552f7358SJed Brown   PetscFunctionBegin;
259552f7358SJed Brown   PetscValidIntPointer(dof, 2);
260552f7358SJed Brown   *dof = ctx->dof;
261552f7358SJed Brown   PetscFunctionReturn(0);
262552f7358SJed Brown }
263552f7358SJed Brown 
2644267b1a3SMatthew G. Knepley /*@C
2654267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2664267b1a3SMatthew G. Knepley 
2674267b1a3SMatthew G. Knepley   Not collective
2684267b1a3SMatthew G. Knepley 
2694267b1a3SMatthew G. Knepley   Input Parameters:
2704267b1a3SMatthew G. Knepley + ctx    - the context
2714267b1a3SMatthew G. Knepley . n      - the number of points
2724267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2734267b1a3SMatthew G. Knepley 
2744267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
2754267b1a3SMatthew G. Knepley 
2764267b1a3SMatthew G. Knepley   Level: intermediate
2774267b1a3SMatthew G. Knepley 
278db781477SPatrick Sanan .seealso: `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()`
2794267b1a3SMatthew G. Knepley @*/
2809371c9d4SSatish Balay PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) {
281552f7358SJed Brown   PetscFunctionBegin;
28208401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
28328b400f6SJacob Faibussowitsch   PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
284552f7358SJed Brown   ctx->nInput = n;
2851aa26658SKarl Rupp 
2869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points));
2879566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim));
288552f7358SJed Brown   PetscFunctionReturn(0);
289552f7358SJed Brown }
290552f7358SJed Brown 
2914267b1a3SMatthew G. Knepley /*@C
29252aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
2934267b1a3SMatthew G. Knepley 
2944267b1a3SMatthew G. Knepley   Collective on ctx
2954267b1a3SMatthew G. Knepley 
2964267b1a3SMatthew G. Knepley   Input Parameters:
2974267b1a3SMatthew G. Knepley + ctx - the context
2984267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
29952aa1562SMatthew G. Knepley . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
3007deeda50SSatish Balay - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
3014267b1a3SMatthew G. Knepley 
3024267b1a3SMatthew G. Knepley   Level: intermediate
3034267b1a3SMatthew G. Knepley 
304db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
3054267b1a3SMatthew G. Knepley @*/
3069371c9d4SSatish Balay PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain) {
307552f7358SJed Brown   MPI_Comm           comm = ctx->comm;
308552f7358SJed Brown   PetscScalar       *a;
309552f7358SJed Brown   PetscInt           p, q, i;
310552f7358SJed Brown   PetscMPIInt        rank, size;
311552f7358SJed Brown   Vec                pointVec;
3123a93e3b7SToby Isaac   PetscSF            cellSF;
313552f7358SJed Brown   PetscLayout        layout;
314552f7358SJed Brown   PetscReal         *globalPoints;
315cb313848SJed Brown   PetscScalar       *globalPointsScalar;
316552f7358SJed Brown   const PetscInt    *ranges;
317552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3183a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3193a93e3b7SToby Isaac   const PetscInt    *foundPoints;
320552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3213a93e3b7SToby Isaac   PetscInt           n, N, numFound;
322552f7358SJed Brown 
32319436ca2SJed Brown   PetscFunctionBegin;
324064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
3259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
32708401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
32819436ca2SJed Brown   /* Locate points */
32919436ca2SJed Brown   n = ctx->nInput;
330552f7358SJed Brown   if (!redundantPoints) {
3319566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &layout));
3329566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(layout, 1));
3339566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(layout, n));
3349566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(layout));
3359566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(layout, &N));
336552f7358SJed Brown     /* Communicate all points to all processes */
3379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs));
3389566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(layout, &ranges));
339552f7358SJed Brown     for (p = 0; p < size; ++p) {
340552f7358SJed Brown       counts[p] = (ranges[p + 1] - ranges[p]) * ctx->dim;
341552f7358SJed Brown       displs[p] = ranges[p] * ctx->dim;
342552f7358SJed Brown     }
3439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(ctx->points, n * ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm));
344552f7358SJed Brown   } else {
345552f7358SJed Brown     N            = n;
346552f7358SJed Brown     globalPoints = ctx->points;
34738ea73c8SJed Brown     counts = displs = NULL;
34838ea73c8SJed Brown     layout          = NULL;
349552f7358SJed Brown   }
350552f7358SJed Brown #if 0
3519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs));
35219436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
353552f7358SJed Brown #else
354cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
3559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar));
35646b3086cSToby Isaac   for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
357cb313848SJed Brown #else
358cb313848SJed Brown   globalPointsScalar = globalPoints;
359cb313848SJed Brown #endif
3609566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec));
3619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs));
3627b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) { foundProcs[p] = size; }
3633a93e3b7SToby Isaac   cellSF = NULL;
3649566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF));
3659566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells));
366552f7358SJed Brown #endif
3673a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3683a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
369552f7358SJed Brown   }
370552f7358SJed Brown   /* Let the lowest rank process own each point */
3711c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm));
372552f7358SJed Brown   ctx->n = 0;
373552f7358SJed Brown   for (p = 0; p < N; ++p) {
37452aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
3759371c9d4SSatish 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),
3769371c9d4SSatish Balay                  (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0.0));
377f7d195e4SLawrence Mitchell       if (rank == 0) ++ctx->n;
37852aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
379552f7358SJed Brown   }
380552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
3819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->n, &ctx->cells));
3829566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &ctx->coords));
3839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE));
3849566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(ctx->coords, ctx->dim));
3859566063dSJacob Faibussowitsch   PetscCall(VecSetType(ctx->coords, VECSTANDARD));
3869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->coords, &a));
387552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
388552f7358SJed Brown     if (globalProcs[p] == rank) {
389552f7358SJed Brown       PetscInt d;
390552f7358SJed Brown 
3911aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d];
392f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
393f317b747SMatthew G. Knepley       ++q;
394552f7358SJed Brown     }
395dd400576SPatrick Sanan     if (globalProcs[p] == size && rank == 0) {
39652aa1562SMatthew G. Knepley       PetscInt d;
39752aa1562SMatthew G. Knepley 
39852aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
39952aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
40052aa1562SMatthew G. Knepley       ++q;
40152aa1562SMatthew G. Knepley     }
402552f7358SJed Brown   }
4039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->coords, &a));
404552f7358SJed Brown #if 0
4059566063dSJacob Faibussowitsch   PetscCall(PetscFree3(foundCells,foundProcs,globalProcs));
406552f7358SJed Brown #else
4079566063dSJacob Faibussowitsch   PetscCall(PetscFree2(foundProcs, globalProcs));
4089566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
4099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pointVec));
410552f7358SJed Brown #endif
4119566063dSJacob Faibussowitsch   if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar));
4129566063dSJacob Faibussowitsch   if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs));
4139566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
414552f7358SJed Brown   PetscFunctionReturn(0);
415552f7358SJed Brown }
416552f7358SJed Brown 
4174267b1a3SMatthew G. Knepley /*@C
4184267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
4194267b1a3SMatthew G. Knepley 
4204267b1a3SMatthew G. Knepley   Collective on ctx
4214267b1a3SMatthew G. Knepley 
4224267b1a3SMatthew G. Knepley   Input Parameter:
4234267b1a3SMatthew G. Knepley . ctx - the context
4244267b1a3SMatthew G. Knepley 
4254267b1a3SMatthew G. Knepley   Output Parameter:
4264267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4274267b1a3SMatthew G. Knepley 
4284267b1a3SMatthew G. Knepley   Note: The local vector entries correspond to interpolation points lying on this process, according to the associated DM. This is a borrowed vector that the user should not destroy.
4294267b1a3SMatthew G. Knepley 
4304267b1a3SMatthew G. Knepley   Level: intermediate
4314267b1a3SMatthew G. Knepley 
432db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4334267b1a3SMatthew G. Knepley @*/
4349371c9d4SSatish Balay PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) {
435552f7358SJed Brown   PetscFunctionBegin;
436552f7358SJed Brown   PetscValidPointer(coordinates, 2);
43728b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
438552f7358SJed Brown   *coordinates = ctx->coords;
439552f7358SJed Brown   PetscFunctionReturn(0);
440552f7358SJed Brown }
441552f7358SJed Brown 
4424267b1a3SMatthew G. Knepley /*@C
4434267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
4444267b1a3SMatthew G. Knepley 
4454267b1a3SMatthew G. Knepley   Collective on ctx
4464267b1a3SMatthew G. Knepley 
4474267b1a3SMatthew G. Knepley   Input Parameter:
4484267b1a3SMatthew G. Knepley . ctx - the context
4494267b1a3SMatthew G. Knepley 
4504267b1a3SMatthew G. Knepley   Output Parameter:
4514267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4524267b1a3SMatthew G. Knepley 
4534267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
4544267b1a3SMatthew G. Knepley 
4554267b1a3SMatthew G. Knepley   Level: intermediate
4564267b1a3SMatthew G. Knepley 
457db781477SPatrick Sanan .seealso: `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4584267b1a3SMatthew G. Knepley @*/
4599371c9d4SSatish Balay PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) {
460552f7358SJed Brown   PetscFunctionBegin;
461552f7358SJed Brown   PetscValidPointer(v, 2);
46228b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4639566063dSJacob Faibussowitsch   PetscCall(VecCreate(ctx->comm, v));
4649566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE));
4659566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*v, ctx->dof));
4669566063dSJacob Faibussowitsch   PetscCall(VecSetType(*v, VECSTANDARD));
467552f7358SJed Brown   PetscFunctionReturn(0);
468552f7358SJed Brown }
469552f7358SJed Brown 
4704267b1a3SMatthew G. Knepley /*@C
4714267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
4724267b1a3SMatthew G. Knepley 
4734267b1a3SMatthew G. Knepley   Collective on ctx
4744267b1a3SMatthew G. Knepley 
4754267b1a3SMatthew G. Knepley   Input Parameters:
4764267b1a3SMatthew G. Knepley + ctx - the context
4774267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
4784267b1a3SMatthew G. Knepley 
4794267b1a3SMatthew G. Knepley   Level: intermediate
4804267b1a3SMatthew G. Knepley 
481db781477SPatrick Sanan .seealso: `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4824267b1a3SMatthew G. Knepley @*/
4839371c9d4SSatish Balay PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) {
484552f7358SJed Brown   PetscFunctionBegin;
485552f7358SJed Brown   PetscValidPointer(v, 2);
48628b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(v));
488552f7358SJed Brown   PetscFunctionReturn(0);
489552f7358SJed Brown }
490552f7358SJed Brown 
4919371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
492cfe5744fSMatthew G. Knepley   PetscReal          v0, J, invJ, detJ;
493cfe5744fSMatthew G. Knepley   const PetscInt     dof = ctx->dof;
494cfe5744fSMatthew G. Knepley   const PetscScalar *coords;
495cfe5744fSMatthew G. Knepley   PetscScalar       *a;
496cfe5744fSMatthew G. Knepley   PetscInt           p;
497cfe5744fSMatthew G. Knepley 
498cfe5744fSMatthew G. Knepley   PetscFunctionBegin;
4999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
501cfe5744fSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
502cfe5744fSMatthew G. Knepley     PetscInt     c = ctx->cells[p];
503cfe5744fSMatthew G. Knepley     PetscScalar *x = NULL;
504cfe5744fSMatthew G. Knepley     PetscReal    xir[1];
505cfe5744fSMatthew G. Knepley     PetscInt     xSize, comp;
506cfe5744fSMatthew G. Knepley 
5079566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ));
508cfe5744fSMatthew G. Knepley     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
509cfe5744fSMatthew G. Knepley     xir[0] = invJ * PetscRealPart(coords[p] - v0);
5109566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
511cfe5744fSMatthew G. Knepley     if (2 * dof == xSize) {
512cfe5744fSMatthew 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];
513cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
514cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
515cfe5744fSMatthew 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);
5169566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
517cfe5744fSMatthew G. Knepley   }
5189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
520cfe5744fSMatthew G. Knepley   PetscFunctionReturn(0);
521cfe5744fSMatthew G. Knepley }
522cfe5744fSMatthew G. Knepley 
5239371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
524552f7358SJed Brown   PetscReal         *v0, *J, *invJ, detJ;
52556044e6dSMatthew G. Knepley   const PetscScalar *coords;
52656044e6dSMatthew G. Knepley   PetscScalar       *a;
527552f7358SJed Brown   PetscInt           p;
528552f7358SJed Brown 
529552f7358SJed Brown   PetscFunctionBegin;
5309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
533552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
534552f7358SJed Brown     PetscInt     c = ctx->cells[p];
535a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
536552f7358SJed Brown     PetscReal    xi[4];
537552f7358SJed Brown     PetscInt     d, f, comp;
538552f7358SJed Brown 
5399566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
54063a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5419566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5421aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
5431aa26658SKarl Rupp 
544552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
545552f7358SJed Brown       xi[d] = 0.0;
5461aa26658SKarl 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]);
5471aa26658SKarl 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];
548552f7358SJed Brown     }
5499566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
550552f7358SJed Brown   }
5519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5539566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
554552f7358SJed Brown   PetscFunctionReturn(0);
555552f7358SJed Brown }
556552f7358SJed Brown 
5579371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
5587a1931ceSMatthew G. Knepley   PetscReal         *v0, *J, *invJ, detJ;
55956044e6dSMatthew G. Knepley   const PetscScalar *coords;
56056044e6dSMatthew G. Knepley   PetscScalar       *a;
5617a1931ceSMatthew G. Knepley   PetscInt           p;
5627a1931ceSMatthew G. Knepley 
5637a1931ceSMatthew G. Knepley   PetscFunctionBegin;
5649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ));
5659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5677a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5687a1931ceSMatthew G. Knepley     PetscInt       c        = ctx->cells[p];
5697a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5702584bbe8SMatthew G. Knepley     PetscScalar   *x        = NULL;
5717a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5727a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5737a1931ceSMatthew G. Knepley 
5749566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
57563a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5769566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5777a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
5787a1931ceSMatthew G. Knepley 
5797a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5807a1931ceSMatthew G. Knepley       xi[d] = 0.0;
5817a1931ceSMatthew 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]);
5827a1931ceSMatthew 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];
5837a1931ceSMatthew G. Knepley     }
5849566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
5857a1931ceSMatthew G. Knepley   }
5869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5889566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
5897a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
5907a1931ceSMatthew G. Knepley }
5917a1931ceSMatthew G. Knepley 
5929371c9d4SSatish Balay static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) {
593552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
594552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
595552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
596552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
597552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
598552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
599552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
600552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
601552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
602552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
603552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
604552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
605552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
606552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
607552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
60856044e6dSMatthew G. Knepley   const PetscScalar *ref;
60956044e6dSMatthew G. Knepley   PetscScalar       *real;
610552f7358SJed Brown 
611552f7358SJed Brown   PetscFunctionBegin;
6129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
6139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
614552f7358SJed Brown   {
615552f7358SJed Brown     const PetscScalar p0 = ref[0];
616552f7358SJed Brown     const PetscScalar p1 = ref[1];
617552f7358SJed Brown 
618552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
619552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
620552f7358SJed Brown   }
6219566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(28));
6229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
624552f7358SJed Brown   PetscFunctionReturn(0);
625552f7358SJed Brown }
626552f7358SJed Brown 
627af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
6289371c9d4SSatish Balay static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) {
629552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
630552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
631552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
632552f7358SJed Brown   const PetscScalar  x1       = vertices[2];
633552f7358SJed Brown   const PetscScalar  y1       = vertices[3];
634552f7358SJed Brown   const PetscScalar  x2       = vertices[4];
635552f7358SJed Brown   const PetscScalar  y2       = vertices[5];
636552f7358SJed Brown   const PetscScalar  x3       = vertices[6];
637552f7358SJed Brown   const PetscScalar  y3       = vertices[7];
638552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
639552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
64056044e6dSMatthew G. Knepley   const PetscScalar *ref;
641552f7358SJed Brown 
642552f7358SJed Brown   PetscFunctionBegin;
6439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
644552f7358SJed Brown   {
645552f7358SJed Brown     const PetscScalar x       = ref[0];
646552f7358SJed Brown     const PetscScalar y       = ref[1];
647552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
648da80777bSKarl Rupp     PetscScalar       values[4];
649da80777bSKarl Rupp 
6509371c9d4SSatish Balay     values[0] = (x1 - x0 + f_01 * y) * 0.5;
6519371c9d4SSatish Balay     values[1] = (x3 - x0 + f_01 * x) * 0.5;
6529371c9d4SSatish Balay     values[2] = (y1 - y0 + g_01 * y) * 0.5;
6539371c9d4SSatish Balay     values[3] = (y3 - y0 + g_01 * x) * 0.5;
6549566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES));
655552f7358SJed Brown   }
6569566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(30));
6579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
6589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
660552f7358SJed Brown   PetscFunctionReturn(0);
661552f7358SJed Brown }
662552f7358SJed Brown 
6639371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
664fafc0619SMatthew G Knepley   DM                 dmCoord;
665de73a395SMatthew G. Knepley   PetscFE            fem = NULL;
666552f7358SJed Brown   SNES               snes;
667552f7358SJed Brown   KSP                ksp;
668552f7358SJed Brown   PC                 pc;
669552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
670552f7358SJed Brown   Mat                J;
671716009bfSMatthew G. Knepley   PetscTabulation    T = NULL;
67256044e6dSMatthew G. Knepley   const PetscScalar *coords;
67356044e6dSMatthew G. Knepley   PetscScalar       *a;
674716009bfSMatthew G. Knepley   PetscReal          xir[2] = {0., 0.};
675de73a395SMatthew G. Knepley   PetscInt           Nf, p;
6765509d985SMatthew G. Knepley   const PetscInt     dof = ctx->dof;
677552f7358SJed Brown 
678552f7358SJed Brown   PetscFunctionBegin;
6799566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
680716009bfSMatthew G. Knepley   if (Nf) {
681cfe5744fSMatthew G. Knepley     PetscObject  obj;
682cfe5744fSMatthew G. Knepley     PetscClassId id;
683cfe5744fSMatthew G. Knepley 
6849566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, 0, NULL, &obj));
6859566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
6869371c9d4SSatish Balay     if (id == PETSCFE_CLASSID) {
6879371c9d4SSatish Balay       fem = (PetscFE)obj;
6889371c9d4SSatish Balay       PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T));
6899371c9d4SSatish Balay     }
690716009bfSMatthew G. Knepley   }
6919566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
6929566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
6939566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
6949566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_"));
6959566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
6969566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 2, 2));
6979566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
6989566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
6999566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
7009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
7019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 2, 2, 2, 2));
7029566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
7039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
7049566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL));
7059566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL));
7069566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7079566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
7089566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
7099566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
710552f7358SJed Brown 
7119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
7129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
713552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
714a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
715552f7358SJed Brown     PetscScalar *xi;
716552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
717552f7358SJed Brown 
718552f7358SJed Brown     /* Can make this do all points at once */
7199566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
72063a3b9bcSJacob Faibussowitsch     PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2);
7219566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
7229566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
7239566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
7249566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
725552f7358SJed Brown     xi[0] = coords[p * ctx->dim + 0];
726552f7358SJed Brown     xi[1] = coords[p * ctx->dim + 1];
7279566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
7289566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
7299566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
730cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
731cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
732cfe5744fSMatthew G. Knepley     if (4 * dof == xSize) {
7339371c9d4SSatish 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];
734cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
735cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp];
736cfe5744fSMatthew G. Knepley     } else {
7375509d985SMatthew G. Knepley       PetscInt d;
7381aa26658SKarl Rupp 
7392c71b3e2SJacob Faibussowitsch       PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE");
7409371c9d4SSatish Balay       xir[0] = 2.0 * xir[0] - 1.0;
7419371c9d4SSatish Balay       xir[1] = 2.0 * xir[1] - 1.0;
7429566063dSJacob Faibussowitsch       PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T));
7435509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7445509d985SMatthew G. Knepley         a[p * dof + comp] = 0.0;
7459371c9d4SSatish Balay         for (d = 0; d < xSize / dof; ++d) { a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp]; }
7465509d985SMatthew G. Knepley       }
7475509d985SMatthew G. Knepley     }
7489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
7499566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
7509566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
751552f7358SJed Brown   }
7529566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&T));
7539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
7549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
755552f7358SJed Brown 
7569566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
7599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
7609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
761552f7358SJed Brown   PetscFunctionReturn(0);
762552f7358SJed Brown }
763552f7358SJed Brown 
7649371c9d4SSatish Balay static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) {
765552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
766552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
767552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
768552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
7697a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
7707a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
7717a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
772552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
773552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
774552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
7757a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
7767a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
7777a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
778552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
779552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
780552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
781552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
782552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
783552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
784552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
785552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
786552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
787552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
788552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
789552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
790552f7358SJed Brown   const PetscScalar  f_1      = x1 - x0;
791552f7358SJed Brown   const PetscScalar  g_1      = y1 - y0;
792552f7358SJed Brown   const PetscScalar  h_1      = z1 - z0;
793552f7358SJed Brown   const PetscScalar  f_3      = x3 - x0;
794552f7358SJed Brown   const PetscScalar  g_3      = y3 - y0;
795552f7358SJed Brown   const PetscScalar  h_3      = z3 - z0;
796552f7358SJed Brown   const PetscScalar  f_4      = x4 - x0;
797552f7358SJed Brown   const PetscScalar  g_4      = y4 - y0;
798552f7358SJed Brown   const PetscScalar  h_4      = z4 - z0;
799552f7358SJed Brown   const PetscScalar  f_01     = x2 - x1 - x3 + x0;
800552f7358SJed Brown   const PetscScalar  g_01     = y2 - y1 - y3 + y0;
801552f7358SJed Brown   const PetscScalar  h_01     = z2 - z1 - z3 + z0;
802552f7358SJed Brown   const PetscScalar  f_12     = x7 - x3 - x4 + x0;
803552f7358SJed Brown   const PetscScalar  g_12     = y7 - y3 - y4 + y0;
804552f7358SJed Brown   const PetscScalar  h_12     = z7 - z3 - z4 + z0;
805552f7358SJed Brown   const PetscScalar  f_02     = x5 - x1 - x4 + x0;
806552f7358SJed Brown   const PetscScalar  g_02     = y5 - y1 - y4 + y0;
807552f7358SJed Brown   const PetscScalar  h_02     = z5 - z1 - z4 + z0;
808552f7358SJed Brown   const PetscScalar  f_012    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
809552f7358SJed Brown   const PetscScalar  g_012    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
810552f7358SJed Brown   const PetscScalar  h_012    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
81156044e6dSMatthew G. Knepley   const PetscScalar *ref;
81256044e6dSMatthew G. Knepley   PetscScalar       *real;
813552f7358SJed Brown 
814552f7358SJed Brown   PetscFunctionBegin;
8159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
8169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
817552f7358SJed Brown   {
818552f7358SJed Brown     const PetscScalar p0 = ref[0];
819552f7358SJed Brown     const PetscScalar p1 = ref[1];
820552f7358SJed Brown     const PetscScalar p2 = ref[2];
821552f7358SJed Brown 
822552f7358SJed 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;
823552f7358SJed 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;
824552f7358SJed 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;
825552f7358SJed Brown   }
8269566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(114));
8279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
8289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
829552f7358SJed Brown   PetscFunctionReturn(0);
830552f7358SJed Brown }
831552f7358SJed Brown 
8329371c9d4SSatish Balay static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) {
833552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar *)ctx;
834552f7358SJed Brown   const PetscScalar  x0       = vertices[0];
835552f7358SJed Brown   const PetscScalar  y0       = vertices[1];
836552f7358SJed Brown   const PetscScalar  z0       = vertices[2];
8377a1931ceSMatthew G. Knepley   const PetscScalar  x1       = vertices[9];
8387a1931ceSMatthew G. Knepley   const PetscScalar  y1       = vertices[10];
8397a1931ceSMatthew G. Knepley   const PetscScalar  z1       = vertices[11];
840552f7358SJed Brown   const PetscScalar  x2       = vertices[6];
841552f7358SJed Brown   const PetscScalar  y2       = vertices[7];
842552f7358SJed Brown   const PetscScalar  z2       = vertices[8];
8437a1931ceSMatthew G. Knepley   const PetscScalar  x3       = vertices[3];
8447a1931ceSMatthew G. Knepley   const PetscScalar  y3       = vertices[4];
8457a1931ceSMatthew G. Knepley   const PetscScalar  z3       = vertices[5];
846552f7358SJed Brown   const PetscScalar  x4       = vertices[12];
847552f7358SJed Brown   const PetscScalar  y4       = vertices[13];
848552f7358SJed Brown   const PetscScalar  z4       = vertices[14];
849552f7358SJed Brown   const PetscScalar  x5       = vertices[15];
850552f7358SJed Brown   const PetscScalar  y5       = vertices[16];
851552f7358SJed Brown   const PetscScalar  z5       = vertices[17];
852552f7358SJed Brown   const PetscScalar  x6       = vertices[18];
853552f7358SJed Brown   const PetscScalar  y6       = vertices[19];
854552f7358SJed Brown   const PetscScalar  z6       = vertices[20];
855552f7358SJed Brown   const PetscScalar  x7       = vertices[21];
856552f7358SJed Brown   const PetscScalar  y7       = vertices[22];
857552f7358SJed Brown   const PetscScalar  z7       = vertices[23];
858552f7358SJed Brown   const PetscScalar  f_xy     = x2 - x1 - x3 + x0;
859552f7358SJed Brown   const PetscScalar  g_xy     = y2 - y1 - y3 + y0;
860552f7358SJed Brown   const PetscScalar  h_xy     = z2 - z1 - z3 + z0;
861552f7358SJed Brown   const PetscScalar  f_yz     = x7 - x3 - x4 + x0;
862552f7358SJed Brown   const PetscScalar  g_yz     = y7 - y3 - y4 + y0;
863552f7358SJed Brown   const PetscScalar  h_yz     = z7 - z3 - z4 + z0;
864552f7358SJed Brown   const PetscScalar  f_xz     = x5 - x1 - x4 + x0;
865552f7358SJed Brown   const PetscScalar  g_xz     = y5 - y1 - y4 + y0;
866552f7358SJed Brown   const PetscScalar  h_xz     = z5 - z1 - z4 + z0;
867552f7358SJed Brown   const PetscScalar  f_xyz    = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
868552f7358SJed Brown   const PetscScalar  g_xyz    = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
869552f7358SJed Brown   const PetscScalar  h_xyz    = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
87056044e6dSMatthew G. Knepley   const PetscScalar *ref;
871552f7358SJed Brown 
872552f7358SJed Brown   PetscFunctionBegin;
8739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref, &ref));
874552f7358SJed Brown   {
875552f7358SJed Brown     const PetscScalar x       = ref[0];
876552f7358SJed Brown     const PetscScalar y       = ref[1];
877552f7358SJed Brown     const PetscScalar z       = ref[2];
878552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
879da80777bSKarl Rupp     PetscScalar       values[9];
880da80777bSKarl Rupp 
881da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0;
882da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0;
883da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0;
884da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0;
885da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0;
886da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0;
887da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0;
888da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0;
889da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0;
8901aa26658SKarl Rupp 
8919566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES));
892552f7358SJed Brown   }
8939566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(152));
8949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref, &ref));
8959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
8969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
897552f7358SJed Brown   PetscFunctionReturn(0);
898552f7358SJed Brown }
899552f7358SJed Brown 
9009371c9d4SSatish Balay static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
901fafc0619SMatthew G Knepley   DM                 dmCoord;
902552f7358SJed Brown   SNES               snes;
903552f7358SJed Brown   KSP                ksp;
904552f7358SJed Brown   PC                 pc;
905552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
906552f7358SJed Brown   Mat                J;
90756044e6dSMatthew G. Knepley   const PetscScalar *coords;
90856044e6dSMatthew G. Knepley   PetscScalar       *a;
909552f7358SJed Brown   PetscInt           p;
910552f7358SJed Brown 
911552f7358SJed Brown   PetscFunctionBegin;
9129566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
9139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
9149566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
9159566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_"));
9169566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
9179566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 3, 3));
9189566063dSJacob Faibussowitsch   PetscCall(VecSetType(r, dm->vectype));
9199566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
9209566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
9219566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
9229566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 3, 3, 3, 3));
9239566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
9249566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
9259566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL));
9269566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL));
9279566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
9289566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
9299566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
9309566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
931552f7358SJed Brown 
9329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
9339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
934552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
935a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
936552f7358SJed Brown     PetscScalar *xi;
937cb313848SJed Brown     PetscReal    xir[3];
938552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
939552f7358SJed Brown 
940552f7358SJed Brown     /* Can make this do all points at once */
9419566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
942cfe5744fSMatthew G. Knepley     PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3);
9439566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
944cfe5744fSMatthew 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);
9459566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
9469566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
9479566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
948552f7358SJed Brown     xi[0] = coords[p * ctx->dim + 0];
949552f7358SJed Brown     xi[1] = coords[p * ctx->dim + 1];
950552f7358SJed Brown     xi[2] = coords[p * ctx->dim + 2];
9519566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
9529566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
9539566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
954cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
955cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
956cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
957cfe5744fSMatthew G. Knepley     if (8 * ctx->dof == xSize) {
958552f7358SJed Brown       for (comp = 0; comp < ctx->dof; ++comp) {
9599371c9d4SSatish 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]) +
9609371c9d4SSatish 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];
961552f7358SJed Brown       }
962cfe5744fSMatthew G. Knepley     } else {
963cfe5744fSMatthew G. Knepley       for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp];
964cfe5744fSMatthew G. Knepley     }
9659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
9669566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
9679566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
968552f7358SJed Brown   }
9699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
9709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
971552f7358SJed Brown 
9729566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
9739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
9749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
9759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
9769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
977552f7358SJed Brown   PetscFunctionReturn(0);
978552f7358SJed Brown }
979552f7358SJed Brown 
9804267b1a3SMatthew G. Knepley /*@C
9814267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9824267b1a3SMatthew G. Knepley 
983552f7358SJed Brown   Input Parameters:
984552f7358SJed Brown + ctx - The DMInterpolationInfo context
985552f7358SJed Brown . dm  - The DM
986552f7358SJed Brown - x   - The local vector containing the field to be interpolated
987552f7358SJed Brown 
988552f7358SJed Brown   Output Parameters:
989552f7358SJed Brown . v   - The vector containing the interpolated values
9904267b1a3SMatthew G. Knepley 
9914267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
9924267b1a3SMatthew G. Knepley 
9934267b1a3SMatthew G. Knepley   Level: beginner
9944267b1a3SMatthew G. Knepley 
995db781477SPatrick Sanan .seealso: `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
9964267b1a3SMatthew G. Knepley @*/
9979371c9d4SSatish Balay PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) {
99866a0a883SMatthew G. Knepley   PetscDS   ds;
99966a0a883SMatthew G. Knepley   PetscInt  n, p, Nf, field;
100066a0a883SMatthew G. Knepley   PetscBool useDS = PETSC_FALSE;
1001552f7358SJed Brown 
1002552f7358SJed Brown   PetscFunctionBegin;
1003552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1004552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1005552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
10069566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
100763a3b9bcSJacob 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);
100866a0a883SMatthew G. Knepley   if (!n) PetscFunctionReturn(0);
10099566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1010680d18d5SMatthew G. Knepley   if (ds) {
101166a0a883SMatthew G. Knepley     useDS = PETSC_TRUE;
10129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &Nf));
1013680d18d5SMatthew G. Knepley     for (field = 0; field < Nf; ++field) {
101466a0a883SMatthew G. Knepley       PetscObject  obj;
1015680d18d5SMatthew G. Knepley       PetscClassId id;
1016680d18d5SMatthew G. Knepley 
10179566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, field, &obj));
10189566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
10199371c9d4SSatish Balay       if (id != PETSCFE_CLASSID) {
10209371c9d4SSatish Balay         useDS = PETSC_FALSE;
10219371c9d4SSatish Balay         break;
10229371c9d4SSatish Balay       }
102366a0a883SMatthew G. Knepley     }
102466a0a883SMatthew G. Knepley   }
102566a0a883SMatthew G. Knepley   if (useDS) {
102666a0a883SMatthew G. Knepley     const PetscScalar *coords;
102766a0a883SMatthew G. Knepley     PetscScalar       *interpolant;
102866a0a883SMatthew G. Knepley     PetscInt           cdim, d;
102966a0a883SMatthew G. Knepley 
10309566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm, &cdim));
10319566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ctx->coords, &coords));
10329566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(v, &interpolant));
1033680d18d5SMatthew G. Knepley     for (p = 0; p < ctx->n; ++p) {
103466a0a883SMatthew G. Knepley       PetscReal    pcoords[3], xi[3];
1035680d18d5SMatthew G. Knepley       PetscScalar *xa   = NULL;
103666a0a883SMatthew G. Knepley       PetscInt     coff = 0, foff = 0, clSize;
1037680d18d5SMatthew G. Knepley 
103852aa1562SMatthew G. Knepley       if (ctx->cells[p] < 0) continue;
1039680d18d5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]);
10409566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi));
10419566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
104266a0a883SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
104366a0a883SMatthew G. Knepley         PetscTabulation T;
104466a0a883SMatthew G. Knepley         PetscFE         fe;
104566a0a883SMatthew G. Knepley 
10469566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
10479566063dSJacob Faibussowitsch         PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T));
1048680d18d5SMatthew G. Knepley         {
1049680d18d5SMatthew G. Knepley           const PetscReal *basis = T->T[0];
1050680d18d5SMatthew G. Knepley           const PetscInt   Nb    = T->Nb;
1051680d18d5SMatthew G. Knepley           const PetscInt   Nc    = T->Nc;
105266a0a883SMatthew G. Knepley           PetscInt         f, fc;
105366a0a883SMatthew G. Knepley 
1054680d18d5SMatthew G. Knepley           for (fc = 0; fc < Nc; ++fc) {
105566a0a883SMatthew G. Knepley             interpolant[p * ctx->dof + coff + fc] = 0.0;
10569371c9d4SSatish Balay             for (f = 0; f < Nb; ++f) { interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc]; }
1057680d18d5SMatthew G. Knepley           }
105866a0a883SMatthew G. Knepley           coff += Nc;
105966a0a883SMatthew G. Knepley           foff += Nb;
1060680d18d5SMatthew G. Knepley         }
10619566063dSJacob Faibussowitsch         PetscCall(PetscTabulationDestroy(&T));
1062680d18d5SMatthew G. Knepley       }
10639566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
106463a3b9bcSJacob Faibussowitsch       PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof);
106563a3b9bcSJacob Faibussowitsch       PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize);
106666a0a883SMatthew G. Knepley     }
10679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
10689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(v, &interpolant));
106966a0a883SMatthew G. Knepley   } else {
107066a0a883SMatthew G. Knepley     DMPolytopeType ct;
107166a0a883SMatthew G. Knepley 
1072680d18d5SMatthew G. Knepley     /* TODO Check each cell individually */
10739566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, ctx->cells[0], &ct));
1074680d18d5SMatthew G. Knepley     switch (ct) {
10759566063dSJacob Faibussowitsch     case DM_POLYTOPE_SEGMENT: PetscCall(DMInterpolate_Segment_Private(ctx, dm, x, v)); break;
10769566063dSJacob Faibussowitsch     case DM_POLYTOPE_TRIANGLE: PetscCall(DMInterpolate_Triangle_Private(ctx, dm, x, v)); break;
10779566063dSJacob Faibussowitsch     case DM_POLYTOPE_QUADRILATERAL: PetscCall(DMInterpolate_Quad_Private(ctx, dm, x, v)); break;
10789566063dSJacob Faibussowitsch     case DM_POLYTOPE_TETRAHEDRON: PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, x, v)); break;
10799566063dSJacob Faibussowitsch     case DM_POLYTOPE_HEXAHEDRON: PetscCall(DMInterpolate_Hex_Private(ctx, dm, x, v)); break;
1080cfe5744fSMatthew G. Knepley     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]);
1081680d18d5SMatthew G. Knepley     }
1082552f7358SJed Brown   }
1083552f7358SJed Brown   PetscFunctionReturn(0);
1084552f7358SJed Brown }
1085552f7358SJed Brown 
10864267b1a3SMatthew G. Knepley /*@C
10874267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
10884267b1a3SMatthew G. Knepley 
10894267b1a3SMatthew G. Knepley   Collective on ctx
10904267b1a3SMatthew G. Knepley 
10914267b1a3SMatthew G. Knepley   Input Parameter:
10924267b1a3SMatthew G. Knepley . ctx - the context
10934267b1a3SMatthew G. Knepley 
10944267b1a3SMatthew G. Knepley   Level: beginner
10954267b1a3SMatthew G. Knepley 
1096db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
10974267b1a3SMatthew G. Knepley @*/
10989371c9d4SSatish Balay PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) {
1099552f7358SJed Brown   PetscFunctionBegin;
1100064a246eSJacob Faibussowitsch   PetscValidPointer(ctx, 1);
11019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->coords));
11029566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->points));
11039566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->cells));
11049566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
11050298fd71SBarry Smith   *ctx = NULL;
1106552f7358SJed Brown   PetscFunctionReturn(0);
1107552f7358SJed Brown }
1108cc0c4584SMatthew G. Knepley 
1109cc0c4584SMatthew G. Knepley /*@C
1110cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1111cc0c4584SMatthew G. Knepley 
1112cc0c4584SMatthew G. Knepley   Collective on SNES
1113cc0c4584SMatthew G. Knepley 
1114cc0c4584SMatthew G. Knepley   Input Parameters:
1115cc0c4584SMatthew G. Knepley + snes   - the SNES context
1116cc0c4584SMatthew G. Knepley . its    - iteration number
1117cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1118d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
1119cc0c4584SMatthew G. Knepley 
1120cc0c4584SMatthew G. Knepley   Notes:
1121cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1122cc0c4584SMatthew G. Knepley 
1123cc0c4584SMatthew G. Knepley   Level: intermediate
1124cc0c4584SMatthew G. Knepley 
1125db781477SPatrick Sanan .seealso: `SNESMonitorSet()`, `SNESMonitorDefault()`
1126cc0c4584SMatthew G. Knepley @*/
11279371c9d4SSatish Balay PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) {
1128d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1129cc0c4584SMatthew G. Knepley   Vec                res;
1130cc0c4584SMatthew G. Knepley   DM                 dm;
1131cc0c4584SMatthew G. Knepley   PetscSection       s;
1132cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1133cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1134cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1135cc0c4584SMatthew G. Knepley 
1136cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11374d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
11389566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
11399566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
11409566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
11419566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
11429566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
11439566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms));
11449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(res, &r));
1145cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1146cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1147cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1148cc0c4584SMatthew G. Knepley 
11499566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
11509566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
1151cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d]));
1152cc0c4584SMatthew G. Knepley     }
1153cc0c4584SMatthew G. Knepley   }
11549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(res, &r));
11551c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
11569566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
11579566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
115863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm));
1159cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
11609566063dSJacob Faibussowitsch     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
11619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f])));
1162cc0c4584SMatthew G. Knepley   }
11639566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
11649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
11659566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
11669566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lnorms, norms));
1167cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1168cc0c4584SMatthew G. Knepley }
116924cdb843SMatthew G. Knepley 
117024cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
117124cdb843SMatthew G. Knepley 
11729371c9d4SSatish Balay PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS) {
11736528b96dSMatthew G. Knepley   PetscInt depth;
11746528b96dSMatthew G. Knepley 
11756528b96dSMatthew G. Knepley   PetscFunctionBegin;
11769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
11779566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS));
11789566063dSJacob Faibussowitsch   if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS));
11796528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
11806528b96dSMatthew G. Knepley }
11816528b96dSMatthew G. Knepley 
118224cdb843SMatthew G. Knepley /*@
11838db23a53SJed Brown   DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user
118424cdb843SMatthew G. Knepley 
118524cdb843SMatthew G. Knepley   Input Parameters:
118624cdb843SMatthew G. Knepley + dm - The mesh
118724cdb843SMatthew G. Knepley . X  - Local solution
118824cdb843SMatthew G. Knepley - user - The user context
118924cdb843SMatthew G. Knepley 
119024cdb843SMatthew G. Knepley   Output Parameter:
119124cdb843SMatthew G. Knepley . F  - Local output vector
119224cdb843SMatthew G. Knepley 
11938db23a53SJed Brown   Notes:
11948db23a53SJed Brown   The residual is summed into F; the caller is responsible for using VecZeroEntries() or otherwise ensuring that any data in F is intentional.
11958db23a53SJed Brown 
119624cdb843SMatthew G. Knepley   Level: developer
119724cdb843SMatthew G. Knepley 
1198db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianAction()`
119924cdb843SMatthew G. Knepley @*/
12009371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) {
12016da023fcSToby Isaac   DM       plex;
1202083401c6SMatthew G. Knepley   IS       allcellIS;
12036528b96dSMatthew G. Knepley   PetscInt Nds, s;
120424cdb843SMatthew G. Knepley 
120524cdb843SMatthew G. Knepley   PetscFunctionBegin;
12069566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12079566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12089566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
12096528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
12106528b96dSMatthew G. Knepley     PetscDS      ds;
12116528b96dSMatthew G. Knepley     IS           cellIS;
121206ad1575SMatthew G. Knepley     PetscFormKey key;
12136528b96dSMatthew G. Knepley 
12149566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
12156528b96dSMatthew G. Knepley     key.value = 0;
12166528b96dSMatthew G. Knepley     key.field = 0;
121706ad1575SMatthew G. Knepley     key.part  = 0;
12186528b96dSMatthew G. Knepley     if (!key.label) {
12199566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
12206528b96dSMatthew G. Knepley       cellIS = allcellIS;
12216528b96dSMatthew G. Knepley     } else {
12226528b96dSMatthew G. Knepley       IS pointIS;
12236528b96dSMatthew G. Knepley 
12246528b96dSMatthew G. Knepley       key.value = 1;
12259566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
12269566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12279566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
12286528b96dSMatthew G. Knepley     }
12299566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
12309566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
12316528b96dSMatthew G. Knepley   }
12329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
12339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
12346528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
12356528b96dSMatthew G. Knepley }
12366528b96dSMatthew G. Knepley 
12379371c9d4SSatish Balay PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user) {
12386528b96dSMatthew G. Knepley   DM       plex;
12396528b96dSMatthew G. Knepley   IS       allcellIS;
12406528b96dSMatthew G. Knepley   PetscInt Nds, s;
12416528b96dSMatthew G. Knepley 
12426528b96dSMatthew G. Knepley   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12459566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1246083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1247083401c6SMatthew G. Knepley     PetscDS ds;
1248083401c6SMatthew G. Knepley     DMLabel label;
1249083401c6SMatthew G. Knepley     IS      cellIS;
1250083401c6SMatthew G. Knepley 
12519566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
12526528b96dSMatthew G. Knepley     {
125306ad1575SMatthew G. Knepley       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
12546528b96dSMatthew G. Knepley       PetscWeakForm     wf;
12556528b96dSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
125606ad1575SMatthew G. Knepley       PetscFormKey     *reskeys;
12576528b96dSMatthew G. Knepley 
12586528b96dSMatthew G. Knepley       /* Get unique residual keys */
12596528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
12606528b96dSMatthew G. Knepley         PetscInt Nkm;
12619566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm));
12626528b96dSMatthew G. Knepley         Nk += Nkm;
12636528b96dSMatthew G. Knepley       }
12649566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &reskeys));
1265*48a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys));
126663a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
12679566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, reskeys));
12686528b96dSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
12696528b96dSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
12706528b96dSMatthew G. Knepley           ++k;
12716528b96dSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
12726528b96dSMatthew G. Knepley         }
12736528b96dSMatthew G. Knepley       }
12746528b96dSMatthew G. Knepley       Nk = k;
12756528b96dSMatthew G. Knepley 
12769566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
12776528b96dSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
12786528b96dSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
12796528b96dSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
12806528b96dSMatthew G. Knepley 
1281083401c6SMatthew G. Knepley         if (!label) {
12829566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
1283083401c6SMatthew G. Knepley           cellIS = allcellIS;
1284083401c6SMatthew G. Knepley         } else {
1285083401c6SMatthew G. Knepley           IS pointIS;
1286083401c6SMatthew G. Knepley 
12879566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
12889566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12899566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
12904a3e9fdbSToby Isaac         }
12919566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
12929566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
1293083401c6SMatthew G. Knepley       }
12949566063dSJacob Faibussowitsch       PetscCall(PetscFree(reskeys));
12956528b96dSMatthew G. Knepley     }
12966528b96dSMatthew G. Knepley   }
12979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
12989566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
129924cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
130024cdb843SMatthew G. Knepley }
130124cdb843SMatthew G. Knepley 
1302bdd6f66aSToby Isaac /*@
1303bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1304bdd6f66aSToby Isaac 
1305bdd6f66aSToby Isaac   Input Parameters:
1306bdd6f66aSToby Isaac + dm - The mesh
1307bdd6f66aSToby Isaac - user - The user context
1308bdd6f66aSToby Isaac 
1309bdd6f66aSToby Isaac   Output Parameter:
1310bdd6f66aSToby Isaac . X  - Local solution
1311bdd6f66aSToby Isaac 
1312bdd6f66aSToby Isaac   Level: developer
1313bdd6f66aSToby Isaac 
1314db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianAction()`
1315bdd6f66aSToby Isaac @*/
13169371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) {
1317bdd6f66aSToby Isaac   DM plex;
1318bdd6f66aSToby Isaac 
1319bdd6f66aSToby Isaac   PetscFunctionBegin;
13209566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
13219566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL));
13229566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1323bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1324bdd6f66aSToby Isaac }
1325bdd6f66aSToby Isaac 
13267a73cf09SMatthew G. Knepley /*@
13278e3a2eefSMatthew G. Knepley   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
13287a73cf09SMatthew G. Knepley 
13297a73cf09SMatthew G. Knepley   Input Parameters:
13308e3a2eefSMatthew G. Knepley + dm   - The DM
13317a73cf09SMatthew G. Knepley . X    - Local solution vector
13327a73cf09SMatthew G. Knepley . Y    - Local input vector
13337a73cf09SMatthew G. Knepley - user - The user context
13347a73cf09SMatthew G. Knepley 
13357a73cf09SMatthew G. Knepley   Output Parameter:
13368e3a2eefSMatthew G. Knepley . F    - lcoal output vector
13377a73cf09SMatthew G. Knepley 
13387a73cf09SMatthew G. Knepley   Level: developer
13397a73cf09SMatthew G. Knepley 
13408e3a2eefSMatthew G. Knepley   Notes:
13418e3a2eefSMatthew G. Knepley   Users will typically use DMSNESCreateJacobianMF() followed by MatMult() instead of calling this routine directly.
13428e3a2eefSMatthew G. Knepley 
1343db781477SPatrick Sanan .seealso: `DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()`
13447a73cf09SMatthew G. Knepley @*/
13459371c9d4SSatish Balay PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user) {
13468e3a2eefSMatthew G. Knepley   DM       plex;
13478e3a2eefSMatthew G. Knepley   IS       allcellIS;
13488e3a2eefSMatthew G. Knepley   PetscInt Nds, s;
1349a925c78cSMatthew G. Knepley 
1350a925c78cSMatthew G. Knepley   PetscFunctionBegin;
13519566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
13529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
13539566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
13548e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
13558e3a2eefSMatthew G. Knepley     PetscDS ds;
13568e3a2eefSMatthew G. Knepley     DMLabel label;
13578e3a2eefSMatthew G. Knepley     IS      cellIS;
13587a73cf09SMatthew G. Knepley 
13599566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
13608e3a2eefSMatthew G. Knepley     {
136106ad1575SMatthew G. Knepley       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
13628e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
13638e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
136406ad1575SMatthew G. Knepley       PetscFormKey     *jackeys;
13658e3a2eefSMatthew G. Knepley 
13668e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
13678e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
13688e3a2eefSMatthew G. Knepley         PetscInt Nkm;
13699566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm));
13708e3a2eefSMatthew G. Knepley         Nk += Nkm;
13718e3a2eefSMatthew G. Knepley       }
13729566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &jackeys));
1373*48a46eb9SPierre Jolivet       for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys));
137463a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
13759566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, jackeys));
13768e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
13778e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
13788e3a2eefSMatthew G. Knepley           ++k;
13798e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
13808e3a2eefSMatthew G. Knepley         }
13818e3a2eefSMatthew G. Knepley       }
13828e3a2eefSMatthew G. Knepley       Nk = k;
13838e3a2eefSMatthew G. Knepley 
13849566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
13858e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
13868e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
13878e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
13888e3a2eefSMatthew G. Knepley 
13898e3a2eefSMatthew G. Knepley         if (!label) {
13909566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)allcellIS));
13918e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
13927a73cf09SMatthew G. Knepley         } else {
13938e3a2eefSMatthew G. Knepley           IS pointIS;
1394a925c78cSMatthew G. Knepley 
13959566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
13969566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
13979566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
1398a925c78cSMatthew G. Knepley         }
13999566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user));
14009566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
14018e3a2eefSMatthew G. Knepley       }
14029566063dSJacob Faibussowitsch       PetscCall(PetscFree(jackeys));
14038e3a2eefSMatthew G. Knepley     }
14048e3a2eefSMatthew G. Knepley   }
14059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
14069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1407a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1408a925c78cSMatthew G. Knepley }
1409a925c78cSMatthew G. Knepley 
141024cdb843SMatthew G. Knepley /*@
141124cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
141224cdb843SMatthew G. Knepley 
141324cdb843SMatthew G. Knepley   Input Parameters:
141424cdb843SMatthew G. Knepley + dm - The mesh
141524cdb843SMatthew G. Knepley . X  - Local input vector
141624cdb843SMatthew G. Knepley - user - The user context
141724cdb843SMatthew G. Knepley 
141824cdb843SMatthew G. Knepley   Output Parameter:
141924cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
142024cdb843SMatthew G. Knepley 
142124cdb843SMatthew G. Knepley   Note:
142224cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
142324cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
142424cdb843SMatthew G. Knepley 
142524cdb843SMatthew G. Knepley   Level: developer
142624cdb843SMatthew G. Knepley 
1427db781477SPatrick Sanan .seealso: `FormFunctionLocal()`
142824cdb843SMatthew G. Knepley @*/
14299371c9d4SSatish Balay PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user) {
14306da023fcSToby Isaac   DM        plex;
1431083401c6SMatthew G. Knepley   IS        allcellIS;
1432f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
14336528b96dSMatthew G. Knepley   PetscInt  Nds, s;
143424cdb843SMatthew G. Knepley 
143524cdb843SMatthew G. Knepley   PetscFunctionBegin;
14369566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
14379566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
14389566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1439083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1440083401c6SMatthew G. Knepley     PetscDS      ds;
1441083401c6SMatthew G. Knepley     IS           cellIS;
144206ad1575SMatthew G. Knepley     PetscFormKey key;
1443083401c6SMatthew G. Knepley 
14449566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
14456528b96dSMatthew G. Knepley     key.value = 0;
14466528b96dSMatthew G. Knepley     key.field = 0;
144706ad1575SMatthew G. Knepley     key.part  = 0;
14486528b96dSMatthew G. Knepley     if (!key.label) {
14499566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1450083401c6SMatthew G. Knepley       cellIS = allcellIS;
1451083401c6SMatthew G. Knepley     } else {
1452083401c6SMatthew G. Knepley       IS pointIS;
1453083401c6SMatthew G. Knepley 
14546528b96dSMatthew G. Knepley       key.value = 1;
14559566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
14569566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14579566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1458083401c6SMatthew G. Knepley     }
1459083401c6SMatthew G. Knepley     if (!s) {
14609566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
14619566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
14629566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
14639566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
1464083401c6SMatthew G. Knepley     }
14659566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user));
14669566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1467083401c6SMatthew G. Knepley   }
14689566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
14699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
147024cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
147124cdb843SMatthew G. Knepley }
14721878804aSMatthew G. Knepley 
14739371c9d4SSatish Balay struct _DMSNESJacobianMFCtx {
14748e3a2eefSMatthew G. Knepley   DM    dm;
14758e3a2eefSMatthew G. Knepley   Vec   X;
14768e3a2eefSMatthew G. Knepley   void *ctx;
14778e3a2eefSMatthew G. Knepley };
14788e3a2eefSMatthew G. Knepley 
14799371c9d4SSatish Balay static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A) {
14808e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
14818e3a2eefSMatthew G. Knepley 
14828e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
14839566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
14849566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, NULL));
14859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
14869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->X));
14879566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
14888e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
14898e3a2eefSMatthew G. Knepley }
14908e3a2eefSMatthew G. Knepley 
14919371c9d4SSatish Balay static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z) {
14928e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
14938e3a2eefSMatthew G. Knepley 
14948e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
14959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
14969566063dSJacob Faibussowitsch   PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx));
14978e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
14988e3a2eefSMatthew G. Knepley }
14998e3a2eefSMatthew G. Knepley 
15008e3a2eefSMatthew G. Knepley /*@
15018e3a2eefSMatthew G. Knepley   DMSNESCreateJacobianMF - Create a Mat which computes the action of the Jacobian matrix-free
15028e3a2eefSMatthew G. Knepley 
15038e3a2eefSMatthew G. Knepley   Collective on dm
15048e3a2eefSMatthew G. Knepley 
15058e3a2eefSMatthew G. Knepley   Input Parameters:
15068e3a2eefSMatthew G. Knepley + dm   - The DM
15078e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
15088e3a2eefSMatthew G. Knepley - user - A user context, or NULL
15098e3a2eefSMatthew G. Knepley 
15108e3a2eefSMatthew G. Knepley   Output Parameter:
15118e3a2eefSMatthew G. Knepley . J    - The Mat
15128e3a2eefSMatthew G. Knepley 
15138e3a2eefSMatthew G. Knepley   Level: advanced
15148e3a2eefSMatthew G. Knepley 
15158e3a2eefSMatthew G. Knepley   Notes:
15168e3a2eefSMatthew G. Knepley   Vec X is kept in Mat J, so updating X then updates the evaluation point.
15178e3a2eefSMatthew G. Knepley 
1518db781477SPatrick Sanan .seealso: `DMSNESComputeJacobianAction()`
15198e3a2eefSMatthew G. Knepley @*/
15209371c9d4SSatish Balay PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J) {
15218e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15228e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
15238e3a2eefSMatthew G. Knepley 
15248e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
15269566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
15279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
15289566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
15299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, n, n, N, N));
15309566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
15319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)X));
15329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &ctx));
15338e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
15348e3a2eefSMatthew G. Knepley   ctx->X   = X;
15358e3a2eefSMatthew G. Knepley   ctx->ctx = user;
15369566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, ctx));
15379566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private));
15389566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private));
15398e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15408e3a2eefSMatthew G. Knepley }
15418e3a2eefSMatthew G. Knepley 
154238cfc46eSPierre Jolivet /*
154338cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
154438cfc46eSPierre Jolivet 
154538cfc46eSPierre Jolivet    Input Parameters:
154638cfc46eSPierre Jolivet +     X - SNES linearization point
154738cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
154838cfc46eSPierre Jolivet 
154938cfc46eSPierre Jolivet    Output Parameter:
155038cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
155138cfc46eSPierre Jolivet 
155238cfc46eSPierre Jolivet    Level: intermediate
155338cfc46eSPierre Jolivet 
1554db781477SPatrick Sanan .seealso: `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()`
155538cfc46eSPierre Jolivet */
15569371c9d4SSatish Balay static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) {
155738cfc46eSPierre Jolivet   SNES   snes;
155838cfc46eSPierre Jolivet   Mat    pJ;
155938cfc46eSPierre Jolivet   DM     ovldm, origdm;
156038cfc46eSPierre Jolivet   DMSNES sdm;
156138cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM, Vec, void *);
156238cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *);
156338cfc46eSPierre Jolivet   void *bctx, *jctx;
156438cfc46eSPierre Jolivet 
156538cfc46eSPierre Jolivet   PetscFunctionBegin;
15669566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ));
156728b400f6SJacob Faibussowitsch   PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat");
15689566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm));
156928b400f6SJacob Faibussowitsch   PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM");
15709566063dSJacob Faibussowitsch   PetscCall(MatGetDM(pJ, &ovldm));
15719566063dSJacob Faibussowitsch   PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx));
15729566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx));
15739566063dSJacob Faibussowitsch   PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx));
15749566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx));
15759566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes));
157638cfc46eSPierre Jolivet   if (!snes) {
15779566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes));
15789566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, ovldm));
15799566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes));
15809566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)snes));
158138cfc46eSPierre Jolivet   }
15829566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(ovldm, &sdm));
15839566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(X));
1584800f99ffSJeremy L Thompson   {
1585800f99ffSJeremy L Thompson     void *ctx;
1586800f99ffSJeremy L Thompson     PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *);
1587800f99ffSJeremy L Thompson     PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx));
1588800f99ffSJeremy L Thompson     PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx));
1589800f99ffSJeremy L Thompson   }
15909566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(X));
159138cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
159238cfc46eSPierre Jolivet   {
159338cfc46eSPierre Jolivet     Mat locpJ;
159438cfc46eSPierre Jolivet 
15959566063dSJacob Faibussowitsch     PetscCall(MatISGetLocalMat(pJ, &locpJ));
15969566063dSJacob Faibussowitsch     PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN));
159738cfc46eSPierre Jolivet   }
159838cfc46eSPierre Jolivet   PetscFunctionReturn(0);
159938cfc46eSPierre Jolivet }
160038cfc46eSPierre Jolivet 
1601a925c78cSMatthew G. Knepley /*@
16029f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
16039f520fc2SToby Isaac 
16049f520fc2SToby Isaac   Input Parameters:
16059f520fc2SToby Isaac + dm - The DM object
1606dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
16079f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
16089f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
16091a244344SSatish Balay 
16101a244344SSatish Balay   Level: developer
16119f520fc2SToby Isaac @*/
16129371c9d4SSatish Balay PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) {
16139f520fc2SToby Isaac   PetscFunctionBegin;
16149566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, boundaryctx));
16159566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, residualctx));
16169566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, jacobianctx));
16179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex));
16189f520fc2SToby Isaac   PetscFunctionReturn(0);
16199f520fc2SToby Isaac }
16209f520fc2SToby Isaac 
1621f017bcb6SMatthew G. Knepley /*@C
1622f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1623f017bcb6SMatthew G. Knepley 
1624f017bcb6SMatthew G. Knepley   Input Parameters:
1625f017bcb6SMatthew G. Knepley + snes - the SNES object
1626f017bcb6SMatthew G. Knepley . dm   - the DM
1627f2cacb80SMatthew G. Knepley . t    - the time
1628f017bcb6SMatthew G. Knepley . u    - a DM vector
1629f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1630f017bcb6SMatthew G. Knepley 
1631f3c94560SMatthew G. Knepley   Output Parameters:
1632f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1633f3c94560SMatthew G. Knepley 
16347f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
16357f96f943SMatthew G. Knepley 
1636f017bcb6SMatthew G. Knepley   Level: developer
1637f017bcb6SMatthew G. Knepley 
1638db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`, `PetscDSSetExactSolution()`
1639f017bcb6SMatthew G. Knepley @*/
16409371c9d4SSatish Balay PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[]) {
1641f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1642f017bcb6SMatthew G. Knepley   void     **ectxs;
1643f3c94560SMatthew G. Knepley   PetscReal *err;
16447f96f943SMatthew G. Knepley   MPI_Comm   comm;
16457f96f943SMatthew G. Knepley   PetscInt   Nf, f;
16461878804aSMatthew G. Knepley 
16471878804aSMatthew G. Knepley   PetscFunctionBegin;
1648f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1649f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1650064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1651534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
16527f96f943SMatthew G. Knepley 
16539566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
16549566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
16557f96f943SMatthew G. Knepley 
16569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
16579566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
16589566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
16597f96f943SMatthew G. Knepley   {
16607f96f943SMatthew G. Knepley     PetscInt Nds, s;
16617f96f943SMatthew G. Knepley 
16629566063dSJacob Faibussowitsch     PetscCall(DMGetNumDS(dm, &Nds));
1663083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
16647f96f943SMatthew G. Knepley       PetscDS         ds;
1665083401c6SMatthew G. Knepley       DMLabel         label;
1666083401c6SMatthew G. Knepley       IS              fieldIS;
16677f96f943SMatthew G. Knepley       const PetscInt *fields;
16687f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1669083401c6SMatthew G. Knepley 
16709566063dSJacob Faibussowitsch       PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds));
16719566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(ds, &dsNf));
16729566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
1673083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1674083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
16759566063dSJacob Faibussowitsch         PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
1676083401c6SMatthew G. Knepley       }
16779566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
1678083401c6SMatthew G. Knepley     }
1679083401c6SMatthew G. Knepley   }
1680f017bcb6SMatthew G. Knepley   if (Nf > 1) {
16819566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err));
1682f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
16839371c9d4SSatish Balay       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); }
1684b878532fSMatthew G. Knepley     } else if (error) {
1685b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
16861878804aSMatthew G. Knepley     } else {
16879566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: ["));
1688f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
16899566063dSJacob Faibussowitsch         if (f) PetscCall(PetscPrintf(comm, ", "));
16909566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(comm, "%g", (double)err[f]));
16911878804aSMatthew G. Knepley       }
16929566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "]\n"));
1693f017bcb6SMatthew G. Knepley     }
1694f017bcb6SMatthew G. Knepley   } else {
16959566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]));
1696f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
169708401ef6SPierre Jolivet       PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol);
1698b878532fSMatthew G. Knepley     } else if (error) {
1699b878532fSMatthew G. Knepley       error[0] = err[0];
1700f017bcb6SMatthew G. Knepley     } else {
17019566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]));
1702f017bcb6SMatthew G. Knepley     }
1703f017bcb6SMatthew G. Knepley   }
17049566063dSJacob Faibussowitsch   PetscCall(PetscFree3(exacts, ectxs, err));
1705f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1706f017bcb6SMatthew G. Knepley }
1707f017bcb6SMatthew G. Knepley 
1708f017bcb6SMatthew G. Knepley /*@C
1709f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1710f017bcb6SMatthew G. Knepley 
1711f017bcb6SMatthew G. Knepley   Input Parameters:
1712f017bcb6SMatthew G. Knepley + snes - the SNES object
1713f017bcb6SMatthew G. Knepley . dm   - the DM
1714f017bcb6SMatthew G. Knepley . u    - a DM vector
1715f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1716f017bcb6SMatthew G. Knepley 
1717f3c94560SMatthew G. Knepley   Output Parameters:
1718f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1719f3c94560SMatthew G. Knepley 
1720f017bcb6SMatthew G. Knepley   Level: developer
1721f017bcb6SMatthew G. Knepley 
1722db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
1723f017bcb6SMatthew G. Knepley @*/
17249371c9d4SSatish Balay PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) {
1725f017bcb6SMatthew G. Knepley   MPI_Comm  comm;
1726f017bcb6SMatthew G. Knepley   Vec       r;
1727f017bcb6SMatthew G. Knepley   PetscReal res;
1728f017bcb6SMatthew G. Knepley 
1729f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1730f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1731f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1732f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1733534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
17349566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
17359566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
17369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
17379566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, u, r));
17389566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
1739f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
174008401ef6SPierre Jolivet     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
1741b878532fSMatthew G. Knepley   } else if (residual) {
1742b878532fSMatthew G. Knepley     *residual = res;
1743f017bcb6SMatthew G. Knepley   } else {
17449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
17459566063dSJacob Faibussowitsch     PetscCall(VecChop(r, 1.0e-10));
17469566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
17479566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
17489566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
1749f017bcb6SMatthew G. Knepley   }
17509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1751f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1752f017bcb6SMatthew G. Knepley }
1753f017bcb6SMatthew G. Knepley 
1754f017bcb6SMatthew G. Knepley /*@C
1755f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1756f017bcb6SMatthew G. Knepley 
1757f017bcb6SMatthew G. Knepley   Input Parameters:
1758f017bcb6SMatthew G. Knepley + snes - the SNES object
1759f017bcb6SMatthew G. Knepley . dm   - the DM
1760f017bcb6SMatthew G. Knepley . u    - a DM vector
1761f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1762f017bcb6SMatthew G. Knepley 
1763f3c94560SMatthew G. Knepley   Output Parameters:
1764f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1765f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1766f3c94560SMatthew G. Knepley 
1767f017bcb6SMatthew G. Knepley   Level: developer
1768f017bcb6SMatthew G. Knepley 
1769db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
1770f017bcb6SMatthew G. Knepley @*/
17719371c9d4SSatish Balay PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) {
1772f017bcb6SMatthew G. Knepley   MPI_Comm     comm;
1773f017bcb6SMatthew G. Knepley   PetscDS      ds;
1774f017bcb6SMatthew G. Knepley   Mat          J, M;
1775f017bcb6SMatthew G. Knepley   MatNullSpace nullspace;
1776f3c94560SMatthew G. Knepley   PetscReal    slope, intercept;
17777f96f943SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
1778f017bcb6SMatthew G. Knepley 
1779f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1780f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1781f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1782f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1783534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1784064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 6);
17859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
17869566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
1787f017bcb6SMatthew G. Knepley   /* Create and view matrices */
17889566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
17899566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
17909566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
17919566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
1792282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
17939566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
17949566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, M));
17959566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
17969566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
17979566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
17989566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
1799282e7bb4SMatthew G. Knepley   } else {
18009566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, J));
1801282e7bb4SMatthew G. Knepley   }
18029566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
18039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
18049566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
1805f017bcb6SMatthew G. Knepley   /* Check nullspace */
18069566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
1807f017bcb6SMatthew G. Knepley   if (nullspace) {
18081878804aSMatthew G. Knepley     PetscBool isNull;
18099566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
181028b400f6SJacob Faibussowitsch     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
18111878804aSMatthew G. Knepley   }
1812f017bcb6SMatthew G. Knepley   /* Taylor test */
1813f017bcb6SMatthew G. Knepley   {
1814f017bcb6SMatthew G. Knepley     PetscRandom rand;
1815f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1816f3c94560SMatthew G. Knepley     PetscReal   h;
1817f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1818f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1819f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1820f017bcb6SMatthew G. Knepley 
1821f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
18229566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
18239566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
18249566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
18259566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
18269566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
18279566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
1828f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
18299566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
18309566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, u, r));
1831f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
18329371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
18339371c9d4SSatish Balay       ;
18349566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
18359566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
18369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
1837f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
18389566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
1839f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
18409566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, uhat, rhat));
18419566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
18429566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
1843f017bcb6SMatthew G. Knepley 
1844f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1845f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1846f017bcb6SMatthew G. Knepley     }
18479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
18489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
18499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
18509566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
18519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
1852f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1853f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1854f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
1855f017bcb6SMatthew G. Knepley     }
1856f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
18579566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
18589566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
1859f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1860f017bcb6SMatthew G. Knepley     if (tol >= 0) {
18610b121fc5SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
1862b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1863b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1864b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1865f017bcb6SMatthew G. Knepley     } else {
18669566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
18679566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
1868f017bcb6SMatthew G. Knepley     }
1869f017bcb6SMatthew G. Knepley   }
18709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
18711878804aSMatthew G. Knepley   PetscFunctionReturn(0);
18721878804aSMatthew G. Knepley }
18731878804aSMatthew G. Knepley 
18749371c9d4SSatish Balay PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u) {
1875f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
18769566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL));
18779566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL));
18789566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL));
1879f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1880f017bcb6SMatthew G. Knepley }
1881f017bcb6SMatthew G. Knepley 
1882bee9a294SMatthew G. Knepley /*@C
1883bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1884bee9a294SMatthew G. Knepley 
1885bee9a294SMatthew G. Knepley   Input Parameters:
1886bee9a294SMatthew G. Knepley + snes - the SNES object
18877f96f943SMatthew G. Knepley - u    - representative SNES vector
18887f96f943SMatthew G. Knepley 
18897f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1890bee9a294SMatthew G. Knepley 
1891bee9a294SMatthew G. Knepley   Level: developer
1892bee9a294SMatthew G. Knepley @*/
18939371c9d4SSatish Balay PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u) {
18941878804aSMatthew G. Knepley   DM        dm;
18951878804aSMatthew G. Knepley   Vec       sol;
18961878804aSMatthew G. Knepley   PetscBool check;
18971878804aSMatthew G. Knepley 
18981878804aSMatthew G. Knepley   PetscFunctionBegin;
18999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check));
19001878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
19019566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
19029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
19039566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, sol));
19049566063dSJacob Faibussowitsch   PetscCall(DMSNESCheck_Internal(snes, dm, sol));
19059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
1906552f7358SJed Brown   PetscFunctionReturn(0);
1907552f7358SJed Brown }
1908