xref: /petsc/src/snes/utils/dmplexsnes.c (revision 6528b96d098a3a0d8b5eec2f1c205a3c25c0d721)
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>
4a925c78cSMatthew G. Knepley #include <petscblaslapack.h>
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
6af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
7552f7358SJed Brown 
8649ef022SMatthew Knepley static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
9649ef022SMatthew Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
10649ef022SMatthew Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
11649ef022SMatthew Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
12649ef022SMatthew Knepley {
13649ef022SMatthew Knepley   p[0] = u[uOff[1]];
14649ef022SMatthew Knepley }
15649ef022SMatthew Knepley 
16649ef022SMatthew Knepley /*
17649ef022SMatthew 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.
18649ef022SMatthew Knepley 
19649ef022SMatthew Knepley   Collective on SNES
20649ef022SMatthew Knepley 
21649ef022SMatthew Knepley   Input Parameters:
22649ef022SMatthew Knepley + snes      - The SNES
23649ef022SMatthew Knepley . pfield    - The field number for pressure
24649ef022SMatthew Knepley . nullspace - The pressure nullspace
25649ef022SMatthew Knepley . u         - The solution vector
26649ef022SMatthew Knepley - ctx       - An optional user context
27649ef022SMatthew Knepley 
28649ef022SMatthew Knepley   Output Parameter:
29649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
30649ef022SMatthew Knepley 
31649ef022SMatthew Knepley   Notes:
32649ef022SMatthew 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.
33649ef022SMatthew Knepley 
34649ef022SMatthew Knepley   Level: developer
35649ef022SMatthew Knepley 
36649ef022SMatthew Knepley .seealso: SNESConvergedCorrectPressure()
37649ef022SMatthew Knepley */
38649ef022SMatthew Knepley static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
39649ef022SMatthew Knepley {
40649ef022SMatthew Knepley   DM             dm;
41649ef022SMatthew Knepley   PetscDS        ds;
42649ef022SMatthew Knepley   const Vec     *nullvecs;
43649ef022SMatthew Knepley   PetscScalar    pintd, *intc, *intn;
44649ef022SMatthew Knepley   MPI_Comm       comm;
45649ef022SMatthew Knepley   PetscInt       Nf, Nv;
46649ef022SMatthew Knepley   PetscErrorCode ierr;
47649ef022SMatthew Knepley 
48649ef022SMatthew Knepley   PetscFunctionBegin;
49649ef022SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
50649ef022SMatthew Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
51649ef022SMatthew Knepley   if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
52649ef022SMatthew Knepley   if (!nullspace) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
53649ef022SMatthew Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
54649ef022SMatthew Knepley   ierr = PetscDSSetObjective(ds, pfield, pressure_Private);CHKERRQ(ierr);
55649ef022SMatthew Knepley   ierr = MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs);CHKERRQ(ierr);
56649ef022SMatthew Knepley   if (Nv != 1) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %D", Nv);
57649ef022SMatthew Knepley   ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
58649ef022SMatthew Knepley   if (PetscAbsScalar(pintd) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
59649ef022SMatthew Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
60649ef022SMatthew Knepley   ierr = PetscMalloc2(Nf, &intc, Nf, &intn);CHKERRQ(ierr);
61649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx);CHKERRQ(ierr);
62649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
63649ef022SMatthew Knepley   ierr = VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]);CHKERRQ(ierr);
64649ef022SMatthew Knepley #if defined (PETSC_USE_DEBUG)
65649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
66649ef022SMatthew Knepley   if (PetscAbsScalar(intc[pfield]) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g\n", (double) PetscRealPart(intc[pfield]));
67649ef022SMatthew Knepley #endif
68649ef022SMatthew Knepley   ierr = PetscFree2(intc, intn);CHKERRQ(ierr);
69649ef022SMatthew Knepley   PetscFunctionReturn(0);
70649ef022SMatthew Knepley }
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley /*@C
73649ef022SMatthew 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().
74649ef022SMatthew Knepley 
75649ef022SMatthew Knepley    Logically Collective on SNES
76649ef022SMatthew Knepley 
77649ef022SMatthew Knepley    Input Parameters:
78649ef022SMatthew Knepley +  snes - the SNES context
79649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
80649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
81649ef022SMatthew Knepley .  snorm - 2-norm of current step
82649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
83649ef022SMatthew Knepley -  ctx   - Optional user context
84649ef022SMatthew Knepley 
85649ef022SMatthew Knepley    Output Parameter:
86649ef022SMatthew Knepley .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
87649ef022SMatthew Knepley 
88649ef022SMatthew Knepley    Notes:
89649ef022SMatthew Knepley    In order to use this monitor, you must setup 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.
90649ef022SMatthew Knepley 
91649ef022SMatthew Knepley    Level: advanced
92649ef022SMatthew Knepley 
93649ef022SMatthew Knepley .seealso: SNESConvergedDefault(), SNESSetConvergenceTest(), DMSetNullSpaceConstructor()
94649ef022SMatthew Knepley @*/
95649ef022SMatthew Knepley PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
96649ef022SMatthew Knepley {
97649ef022SMatthew Knepley   PetscBool      monitorIntegral = PETSC_FALSE;
98649ef022SMatthew Knepley   PetscErrorCode ierr;
99649ef022SMatthew Knepley 
100649ef022SMatthew Knepley   PetscFunctionBegin;
101649ef022SMatthew Knepley   ierr = SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx);CHKERRQ(ierr);
102649ef022SMatthew Knepley   if (monitorIntegral) {
103649ef022SMatthew Knepley     Mat          J;
104649ef022SMatthew Knepley     Vec          u;
105649ef022SMatthew Knepley     MatNullSpace nullspace;
106649ef022SMatthew Knepley     const Vec   *nullvecs;
107649ef022SMatthew Knepley     PetscScalar  pintd;
108649ef022SMatthew Knepley 
109649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
110649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
111649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
112649ef022SMatthew Knepley     ierr = MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);CHKERRQ(ierr);
113649ef022SMatthew Knepley     ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
114649ef022SMatthew Knepley     ierr = PetscInfo1(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));CHKERRQ(ierr);
115649ef022SMatthew Knepley   }
116649ef022SMatthew Knepley   if (*reason > 0) {
117649ef022SMatthew Knepley     Mat          J;
118649ef022SMatthew Knepley     Vec          u;
119649ef022SMatthew Knepley     MatNullSpace nullspace;
120649ef022SMatthew Knepley     PetscInt     pfield = 1;
121649ef022SMatthew Knepley 
122649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
123649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
124649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
125649ef022SMatthew Knepley     ierr = SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx);CHKERRQ(ierr);
126649ef022SMatthew Knepley   }
127649ef022SMatthew Knepley   PetscFunctionReturn(0);
128649ef022SMatthew Knepley }
129649ef022SMatthew Knepley 
13024cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
13124cdb843SMatthew G. Knepley 
1326da023fcSToby Isaac static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
1336da023fcSToby Isaac {
1346da023fcSToby Isaac   PetscBool      isPlex;
1356da023fcSToby Isaac   PetscErrorCode ierr;
1366da023fcSToby Isaac 
1376da023fcSToby Isaac   PetscFunctionBegin;
1386da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
1396da023fcSToby Isaac   if (isPlex) {
1406da023fcSToby Isaac     *plex = dm;
1416da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
142f7148743SMatthew G. Knepley   } else {
143f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
144f7148743SMatthew G. Knepley     if (!*plex) {
1456da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
146f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
1476da023fcSToby Isaac       if (copy) {
1486da023fcSToby Isaac         PetscInt    i;
1496da023fcSToby Isaac         PetscObject obj;
1506da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
1516da023fcSToby Isaac 
1526da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
1536da023fcSToby Isaac         for (i = 0; i < 3; i++) {
1546da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
1556da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
1566da023fcSToby Isaac         }
1576da023fcSToby Isaac       }
158f7148743SMatthew G. Knepley     } else {
159f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
160f7148743SMatthew G. Knepley     }
1616da023fcSToby Isaac   }
1626da023fcSToby Isaac   PetscFunctionReturn(0);
1636da023fcSToby Isaac }
1646da023fcSToby Isaac 
1654267b1a3SMatthew G. Knepley /*@C
1664267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
1674267b1a3SMatthew G. Knepley 
168d083f849SBarry Smith   Collective
1694267b1a3SMatthew G. Knepley 
1704267b1a3SMatthew G. Knepley   Input Parameter:
1714267b1a3SMatthew G. Knepley . comm - the communicator
1724267b1a3SMatthew G. Knepley 
1734267b1a3SMatthew G. Knepley   Output Parameter:
1744267b1a3SMatthew G. Knepley . ctx - the context
1754267b1a3SMatthew G. Knepley 
1764267b1a3SMatthew G. Knepley   Level: beginner
1774267b1a3SMatthew G. Knepley 
1784267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
1794267b1a3SMatthew G. Knepley @*/
1800adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
1810adebc6cSBarry Smith {
182552f7358SJed Brown   PetscErrorCode ierr;
183552f7358SJed Brown 
184552f7358SJed Brown   PetscFunctionBegin;
185552f7358SJed Brown   PetscValidPointer(ctx, 2);
18695dccacaSBarry Smith   ierr = PetscNew(ctx);CHKERRQ(ierr);
1871aa26658SKarl Rupp 
188552f7358SJed Brown   (*ctx)->comm   = comm;
189552f7358SJed Brown   (*ctx)->dim    = -1;
190552f7358SJed Brown   (*ctx)->nInput = 0;
1910298fd71SBarry Smith   (*ctx)->points = NULL;
1920298fd71SBarry Smith   (*ctx)->cells  = NULL;
193552f7358SJed Brown   (*ctx)->n      = -1;
1940298fd71SBarry Smith   (*ctx)->coords = NULL;
195552f7358SJed Brown   PetscFunctionReturn(0);
196552f7358SJed Brown }
197552f7358SJed Brown 
1984267b1a3SMatthew G. Knepley /*@C
1994267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
2004267b1a3SMatthew G. Knepley 
2014267b1a3SMatthew G. Knepley   Not collective
2024267b1a3SMatthew G. Knepley 
2034267b1a3SMatthew G. Knepley   Input Parameters:
2044267b1a3SMatthew G. Knepley + ctx - the context
2054267b1a3SMatthew G. Knepley - dim - the spatial dimension
2064267b1a3SMatthew G. Knepley 
2074267b1a3SMatthew G. Knepley   Level: intermediate
2084267b1a3SMatthew G. Knepley 
2094267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2104267b1a3SMatthew G. Knepley @*/
2110adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
2120adebc6cSBarry Smith {
213552f7358SJed Brown   PetscFunctionBegin;
214ff1e0c32SBarry Smith   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
215552f7358SJed Brown   ctx->dim = dim;
216552f7358SJed Brown   PetscFunctionReturn(0);
217552f7358SJed Brown }
218552f7358SJed Brown 
2194267b1a3SMatthew G. Knepley /*@C
2204267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2214267b1a3SMatthew G. Knepley 
2224267b1a3SMatthew G. Knepley   Not collective
2234267b1a3SMatthew G. Knepley 
2244267b1a3SMatthew G. Knepley   Input Parameter:
2254267b1a3SMatthew G. Knepley . ctx - the context
2264267b1a3SMatthew G. Knepley 
2274267b1a3SMatthew G. Knepley   Output Parameter:
2284267b1a3SMatthew G. Knepley . dim - the spatial dimension
2294267b1a3SMatthew G. Knepley 
2304267b1a3SMatthew G. Knepley   Level: intermediate
2314267b1a3SMatthew G. Knepley 
2324267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2334267b1a3SMatthew G. Knepley @*/
2340adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
2350adebc6cSBarry Smith {
236552f7358SJed Brown   PetscFunctionBegin;
237552f7358SJed Brown   PetscValidIntPointer(dim, 2);
238552f7358SJed Brown   *dim = ctx->dim;
239552f7358SJed Brown   PetscFunctionReturn(0);
240552f7358SJed Brown }
241552f7358SJed Brown 
2424267b1a3SMatthew G. Knepley /*@C
2434267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets 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 Parameters:
2484267b1a3SMatthew G. Knepley + ctx - the context
2494267b1a3SMatthew G. Knepley - dof - the number of fields
2504267b1a3SMatthew G. Knepley 
2514267b1a3SMatthew G. Knepley   Level: intermediate
2524267b1a3SMatthew G. Knepley 
2534267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2544267b1a3SMatthew G. Knepley @*/
2550adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
2560adebc6cSBarry Smith {
257552f7358SJed Brown   PetscFunctionBegin;
258ff1e0c32SBarry Smith   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
259552f7358SJed Brown   ctx->dof = dof;
260552f7358SJed Brown   PetscFunctionReturn(0);
261552f7358SJed Brown }
262552f7358SJed Brown 
2634267b1a3SMatthew G. Knepley /*@C
2644267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2654267b1a3SMatthew G. Knepley 
2664267b1a3SMatthew G. Knepley   Not collective
2674267b1a3SMatthew G. Knepley 
2684267b1a3SMatthew G. Knepley   Input Parameter:
2694267b1a3SMatthew G. Knepley . ctx - the context
2704267b1a3SMatthew G. Knepley 
2714267b1a3SMatthew G. Knepley   Output Parameter:
2724267b1a3SMatthew G. Knepley . dof - the number of fields
2734267b1a3SMatthew G. Knepley 
2744267b1a3SMatthew G. Knepley   Level: intermediate
2754267b1a3SMatthew G. Knepley 
2764267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2774267b1a3SMatthew G. Knepley @*/
2780adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
2790adebc6cSBarry Smith {
280552f7358SJed Brown   PetscFunctionBegin;
281552f7358SJed Brown   PetscValidIntPointer(dof, 2);
282552f7358SJed Brown   *dof = ctx->dof;
283552f7358SJed Brown   PetscFunctionReturn(0);
284552f7358SJed Brown }
285552f7358SJed Brown 
2864267b1a3SMatthew G. Knepley /*@C
2874267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2884267b1a3SMatthew G. Knepley 
2894267b1a3SMatthew G. Knepley   Not collective
2904267b1a3SMatthew G. Knepley 
2914267b1a3SMatthew G. Knepley   Input Parameters:
2924267b1a3SMatthew G. Knepley + ctx    - the context
2934267b1a3SMatthew G. Knepley . n      - the number of points
2944267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2954267b1a3SMatthew G. Knepley 
2964267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
2974267b1a3SMatthew G. Knepley 
2984267b1a3SMatthew G. Knepley   Level: intermediate
2994267b1a3SMatthew G. Knepley 
3004267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
3014267b1a3SMatthew G. Knepley @*/
3020adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
3030adebc6cSBarry Smith {
304552f7358SJed Brown   PetscErrorCode ierr;
305552f7358SJed Brown 
306552f7358SJed Brown   PetscFunctionBegin;
3070adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
3080adebc6cSBarry Smith   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
309552f7358SJed Brown   ctx->nInput = n;
3101aa26658SKarl Rupp 
311785e854fSJed Brown   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
312580bdb30SBarry Smith   ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr);
313552f7358SJed Brown   PetscFunctionReturn(0);
314552f7358SJed Brown }
315552f7358SJed Brown 
3164267b1a3SMatthew G. Knepley /*@C
31752aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
3184267b1a3SMatthew G. Knepley 
3194267b1a3SMatthew G. Knepley   Collective on ctx
3204267b1a3SMatthew G. Knepley 
3214267b1a3SMatthew G. Knepley   Input Parameters:
3224267b1a3SMatthew G. Knepley + ctx - the context
3234267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
32452aa1562SMatthew G. Knepley . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
3257deeda50SSatish Balay - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
3264267b1a3SMatthew G. Knepley 
3274267b1a3SMatthew G. Knepley   Level: intermediate
3284267b1a3SMatthew G. Knepley 
3294267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3304267b1a3SMatthew G. Knepley @*/
33152aa1562SMatthew G. Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
3320adebc6cSBarry Smith {
333552f7358SJed Brown   MPI_Comm          comm = ctx->comm;
334552f7358SJed Brown   PetscScalar       *a;
335552f7358SJed Brown   PetscInt          p, q, i;
336552f7358SJed Brown   PetscMPIInt       rank, size;
337552f7358SJed Brown   PetscErrorCode    ierr;
338552f7358SJed Brown   Vec               pointVec;
3393a93e3b7SToby Isaac   PetscSF           cellSF;
340552f7358SJed Brown   PetscLayout       layout;
341552f7358SJed Brown   PetscReal         *globalPoints;
342cb313848SJed Brown   PetscScalar       *globalPointsScalar;
343552f7358SJed Brown   const PetscInt    *ranges;
344552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3453a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3463a93e3b7SToby Isaac   const PetscInt    *foundPoints;
347552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3483a93e3b7SToby Isaac   PetscInt          n, N, numFound;
349552f7358SJed Brown 
35019436ca2SJed Brown   PetscFunctionBegin;
35119436ca2SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
352ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
353ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
3540adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
35519436ca2SJed Brown   /* Locate points */
35619436ca2SJed Brown   n = ctx->nInput;
357552f7358SJed Brown   if (!redundantPoints) {
358552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
359552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
360552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
361552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
362552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
363552f7358SJed Brown     /* Communicate all points to all processes */
364dcca6d9dSJed Brown     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
365552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
366552f7358SJed Brown     for (p = 0; p < size; ++p) {
367552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
368552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
369552f7358SJed Brown     }
370ffc4695bSBarry Smith     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRMPI(ierr);
371552f7358SJed Brown   } else {
372552f7358SJed Brown     N = n;
373552f7358SJed Brown     globalPoints = ctx->points;
37438ea73c8SJed Brown     counts = displs = NULL;
37538ea73c8SJed Brown     layout = NULL;
376552f7358SJed Brown   }
377552f7358SJed Brown #if 0
378dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
37919436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
380552f7358SJed Brown #else
381cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
38246b3086cSToby Isaac   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
38346b3086cSToby Isaac   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
384cb313848SJed Brown #else
385cb313848SJed Brown   globalPointsScalar = globalPoints;
386cb313848SJed Brown #endif
38704706141SMatthew G Knepley   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
388dcca6d9dSJed Brown   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
3897b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
3903a93e3b7SToby Isaac   cellSF = NULL;
3912d1fa6caSMatthew G. Knepley   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
3923a93e3b7SToby Isaac   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
393552f7358SJed Brown #endif
3943a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3953a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
396552f7358SJed Brown   }
397552f7358SJed Brown   /* Let the lowest rank process own each point */
398b2566f29SBarry Smith   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
399552f7358SJed Brown   ctx->n = 0;
400552f7358SJed Brown   for (p = 0; p < N; ++p) {
40152aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
40252aa1562SMatthew G. Knepley       if (!ignoreOutsideDomain) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %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), (double)(ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0));
40352aa1562SMatthew G. Knepley       else if (!rank) ++ctx->n;
40452aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
405552f7358SJed Brown   }
406552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
407785e854fSJed Brown   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
408552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
409552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
410552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
411c0dedaeaSBarry Smith   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
412552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
413552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
414552f7358SJed Brown     if (globalProcs[p] == rank) {
415552f7358SJed Brown       PetscInt d;
416552f7358SJed Brown 
4171aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
418f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
419f317b747SMatthew G. Knepley       ++q;
420552f7358SJed Brown     }
42152aa1562SMatthew G. Knepley     if (globalProcs[p] == size && !rank) {
42252aa1562SMatthew G. Knepley       PetscInt d;
42352aa1562SMatthew G. Knepley 
42452aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
42552aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
42652aa1562SMatthew G. Knepley       ++q;
42752aa1562SMatthew G. Knepley     }
428552f7358SJed Brown   }
429552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
430552f7358SJed Brown #if 0
431552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
432552f7358SJed Brown #else
433552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
4343a93e3b7SToby Isaac   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
435552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
436552f7358SJed Brown #endif
437cb313848SJed Brown   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
438d343d804SMatthew G. Knepley   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
439552f7358SJed Brown   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
440552f7358SJed Brown   PetscFunctionReturn(0);
441552f7358SJed Brown }
442552f7358SJed Brown 
4434267b1a3SMatthew G. Knepley /*@C
4444267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
4454267b1a3SMatthew G. Knepley 
4464267b1a3SMatthew G. Knepley   Collective on ctx
4474267b1a3SMatthew G. Knepley 
4484267b1a3SMatthew G. Knepley   Input Parameter:
4494267b1a3SMatthew G. Knepley . ctx - the context
4504267b1a3SMatthew G. Knepley 
4514267b1a3SMatthew G. Knepley   Output Parameter:
4524267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4534267b1a3SMatthew G. Knepley 
4544267b1a3SMatthew 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.
4554267b1a3SMatthew G. Knepley 
4564267b1a3SMatthew G. Knepley   Level: intermediate
4574267b1a3SMatthew G. Knepley 
4584267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4594267b1a3SMatthew G. Knepley @*/
4600adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
4610adebc6cSBarry Smith {
462552f7358SJed Brown   PetscFunctionBegin;
463552f7358SJed Brown   PetscValidPointer(coordinates, 2);
4640adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
465552f7358SJed Brown   *coordinates = ctx->coords;
466552f7358SJed Brown   PetscFunctionReturn(0);
467552f7358SJed Brown }
468552f7358SJed Brown 
4694267b1a3SMatthew G. Knepley /*@C
4704267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
4714267b1a3SMatthew G. Knepley 
4724267b1a3SMatthew G. Knepley   Collective on ctx
4734267b1a3SMatthew G. Knepley 
4744267b1a3SMatthew G. Knepley   Input Parameter:
4754267b1a3SMatthew G. Knepley . ctx - the context
4764267b1a3SMatthew G. Knepley 
4774267b1a3SMatthew G. Knepley   Output Parameter:
4784267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4794267b1a3SMatthew G. Knepley 
4804267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
4814267b1a3SMatthew G. Knepley 
4824267b1a3SMatthew G. Knepley   Level: intermediate
4834267b1a3SMatthew G. Knepley 
4844267b1a3SMatthew G. Knepley .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4854267b1a3SMatthew G. Knepley @*/
4860adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
4870adebc6cSBarry Smith {
488552f7358SJed Brown   PetscErrorCode ierr;
489552f7358SJed Brown 
490552f7358SJed Brown   PetscFunctionBegin;
491552f7358SJed Brown   PetscValidPointer(v, 2);
4920adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
493552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
494552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
495552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
496c0dedaeaSBarry Smith   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
497552f7358SJed Brown   PetscFunctionReturn(0);
498552f7358SJed Brown }
499552f7358SJed Brown 
5004267b1a3SMatthew G. Knepley /*@C
5014267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
5024267b1a3SMatthew G. Knepley 
5034267b1a3SMatthew G. Knepley   Collective on ctx
5044267b1a3SMatthew G. Knepley 
5054267b1a3SMatthew G. Knepley   Input Parameters:
5064267b1a3SMatthew G. Knepley + ctx - the context
5074267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
5084267b1a3SMatthew G. Knepley 
5094267b1a3SMatthew G. Knepley   Level: intermediate
5104267b1a3SMatthew G. Knepley 
5114267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
5124267b1a3SMatthew G. Knepley @*/
5130adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
5140adebc6cSBarry Smith {
515552f7358SJed Brown   PetscErrorCode ierr;
516552f7358SJed Brown 
517552f7358SJed Brown   PetscFunctionBegin;
518552f7358SJed Brown   PetscValidPointer(v, 2);
5190adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
520552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
521552f7358SJed Brown   PetscFunctionReturn(0);
522552f7358SJed Brown }
523552f7358SJed Brown 
5247a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
525a6dfd86eSKarl Rupp {
526552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
52756044e6dSMatthew G. Knepley   const PetscScalar *coords;
52856044e6dSMatthew G. Knepley   PetscScalar    *a;
529552f7358SJed Brown   PetscInt       p;
530552f7358SJed Brown   PetscErrorCode ierr;
531552f7358SJed Brown 
532552f7358SJed Brown   PetscFunctionBegin;
533dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
53456044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
535552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
536552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
537552f7358SJed Brown     PetscInt     c = ctx->cells[p];
538a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
539552f7358SJed Brown     PetscReal    xi[4];
540552f7358SJed Brown     PetscInt     d, f, comp;
541552f7358SJed Brown 
5428e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
543ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5440298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5451aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5461aa26658SKarl Rupp 
547552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
548552f7358SJed Brown       xi[d] = 0.0;
5491aa26658SKarl 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]);
5501aa26658SKarl 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];
551552f7358SJed Brown     }
5520298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
553552f7358SJed Brown   }
554552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
55556044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
556552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
557552f7358SJed Brown   PetscFunctionReturn(0);
558552f7358SJed Brown }
559552f7358SJed Brown 
5607a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
5617a1931ceSMatthew G. Knepley {
5627a1931ceSMatthew G. Knepley   PetscReal      *v0, *J, *invJ, detJ;
56356044e6dSMatthew G. Knepley   const PetscScalar *coords;
56456044e6dSMatthew G. Knepley   PetscScalar    *a;
5657a1931ceSMatthew G. Knepley   PetscInt       p;
5667a1931ceSMatthew G. Knepley   PetscErrorCode ierr;
5677a1931ceSMatthew G. Knepley 
5687a1931ceSMatthew G. Knepley   PetscFunctionBegin;
569dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
57056044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5717a1931ceSMatthew G. Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
5727a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5737a1931ceSMatthew G. Knepley     PetscInt       c = ctx->cells[p];
5747a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5752584bbe8SMatthew G. Knepley     PetscScalar   *x = NULL;
5767a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5777a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5787a1931ceSMatthew G. Knepley 
5798e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
580ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5817a1931ceSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5827a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5837a1931ceSMatthew G. Knepley 
5847a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5857a1931ceSMatthew G. Knepley       xi[d] = 0.0;
5867a1931ceSMatthew 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]);
5877a1931ceSMatthew 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];
5887a1931ceSMatthew G. Knepley     }
5897a1931ceSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5907a1931ceSMatthew G. Knepley   }
5917a1931ceSMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
59256044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5937a1931ceSMatthew G. Knepley   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
5947a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
5957a1931ceSMatthew G. Knepley }
5967a1931ceSMatthew G. Knepley 
5975820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
598552f7358SJed Brown {
599552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
600552f7358SJed Brown   const PetscScalar x0        = vertices[0];
601552f7358SJed Brown   const PetscScalar y0        = vertices[1];
602552f7358SJed Brown   const PetscScalar x1        = vertices[2];
603552f7358SJed Brown   const PetscScalar y1        = vertices[3];
604552f7358SJed Brown   const PetscScalar x2        = vertices[4];
605552f7358SJed Brown   const PetscScalar y2        = vertices[5];
606552f7358SJed Brown   const PetscScalar x3        = vertices[6];
607552f7358SJed Brown   const PetscScalar y3        = vertices[7];
608552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
609552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
610552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
611552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
612552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
613552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
61456044e6dSMatthew G. Knepley   const PetscScalar *ref;
61556044e6dSMatthew G. Knepley   PetscScalar       *real;
616552f7358SJed Brown   PetscErrorCode    ierr;
617552f7358SJed Brown 
618552f7358SJed Brown   PetscFunctionBegin;
61956044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
620552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
621552f7358SJed Brown   {
622552f7358SJed Brown     const PetscScalar p0 = ref[0];
623552f7358SJed Brown     const PetscScalar p1 = ref[1];
624552f7358SJed Brown 
625552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
626552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
627552f7358SJed Brown   }
628552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
62956044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
630552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
631552f7358SJed Brown   PetscFunctionReturn(0);
632552f7358SJed Brown }
633552f7358SJed Brown 
634af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
635d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
636552f7358SJed Brown {
637552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
638552f7358SJed Brown   const PetscScalar x0        = vertices[0];
639552f7358SJed Brown   const PetscScalar y0        = vertices[1];
640552f7358SJed Brown   const PetscScalar x1        = vertices[2];
641552f7358SJed Brown   const PetscScalar y1        = vertices[3];
642552f7358SJed Brown   const PetscScalar x2        = vertices[4];
643552f7358SJed Brown   const PetscScalar y2        = vertices[5];
644552f7358SJed Brown   const PetscScalar x3        = vertices[6];
645552f7358SJed Brown   const PetscScalar y3        = vertices[7];
646552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
647552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
64856044e6dSMatthew G. Knepley   const PetscScalar *ref;
649552f7358SJed Brown   PetscErrorCode    ierr;
650552f7358SJed Brown 
651552f7358SJed Brown   PetscFunctionBegin;
65256044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
653552f7358SJed Brown   {
654552f7358SJed Brown     const PetscScalar x       = ref[0];
655552f7358SJed Brown     const PetscScalar y       = ref[1];
656552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
657da80777bSKarl Rupp     PetscScalar       values[4];
658da80777bSKarl Rupp 
659da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
660da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
66194ab13aaSBarry Smith     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
662552f7358SJed Brown   }
663552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
66456044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
66594ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
66694ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
667552f7358SJed Brown   PetscFunctionReturn(0);
668552f7358SJed Brown }
669552f7358SJed Brown 
670a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
671a6dfd86eSKarl Rupp {
672fafc0619SMatthew G Knepley   DM             dmCoord;
673de73a395SMatthew G. Knepley   PetscFE        fem = NULL;
674552f7358SJed Brown   SNES           snes;
675552f7358SJed Brown   KSP            ksp;
676552f7358SJed Brown   PC             pc;
677552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
678552f7358SJed Brown   Mat            J;
679ef0bb6c7SMatthew G. Knepley   PetscTabulation    T;
68056044e6dSMatthew G. Knepley   const PetscScalar *coords;
68156044e6dSMatthew G. Knepley   PetscScalar    *a;
682ef0bb6c7SMatthew G. Knepley   PetscReal      xir[2];
683de73a395SMatthew G. Knepley   PetscInt       Nf, p;
6845509d985SMatthew G. Knepley   const PetscInt dof = ctx->dof;
685552f7358SJed Brown   PetscErrorCode ierr;
686552f7358SJed Brown 
687552f7358SJed Brown   PetscFunctionBegin;
688de73a395SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
68944a7f3ddSMatthew G. Knepley   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
690552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
691fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
692552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
693552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
694552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
695552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
696c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
697552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
698552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
699552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
700552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
701552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
702552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
7030298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
7040298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
705552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
706552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
707552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
708552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
709552f7358SJed Brown 
71056044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
711552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
712ef0bb6c7SMatthew G. Knepley   ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr);
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 */
7190298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
720ff1e0c32SBarry Smith     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
7210298fd71SBarry Smith     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
7220298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
7230298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
724552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
725552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
726552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
727552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
728552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
729552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
730cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
731cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
7325509d985SMatthew G. Knepley     if (4*dof != xSize) {
7335509d985SMatthew G. Knepley       PetscInt d;
7341aa26658SKarl Rupp 
7355509d985SMatthew G. Knepley       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
736ef0bb6c7SMatthew G. Knepley       ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr);
7375509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7385509d985SMatthew G. Knepley         a[p*dof+comp] = 0.0;
7395509d985SMatthew G. Knepley         for (d = 0; d < xSize/dof; ++d) {
740ef0bb6c7SMatthew G. Knepley           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
7415509d985SMatthew G. Knepley         }
7425509d985SMatthew G. Knepley       }
7435509d985SMatthew G. Knepley     } else {
7445509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp)
7455509d985SMatthew G. Knepley         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];
7465509d985SMatthew G. Knepley     }
747552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
7480298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
7490298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
750552f7358SJed Brown   }
751ef0bb6c7SMatthew G. Knepley   ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
752552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
75356044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
754552f7358SJed Brown 
755552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
756552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
757552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
758552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
759552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
760552f7358SJed Brown   PetscFunctionReturn(0);
761552f7358SJed Brown }
762552f7358SJed Brown 
7635820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
764552f7358SJed Brown {
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   PetscErrorCode    ierr;
814552f7358SJed Brown 
815552f7358SJed Brown   PetscFunctionBegin;
81656044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
817552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
818552f7358SJed Brown   {
819552f7358SJed Brown     const PetscScalar p0 = ref[0];
820552f7358SJed Brown     const PetscScalar p1 = ref[1];
821552f7358SJed Brown     const PetscScalar p2 = ref[2];
822552f7358SJed Brown 
823552f7358SJed 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;
824552f7358SJed 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;
825552f7358SJed 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;
826552f7358SJed Brown   }
827552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
82856044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
829552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
830552f7358SJed Brown   PetscFunctionReturn(0);
831552f7358SJed Brown }
832552f7358SJed Brown 
833d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
834552f7358SJed Brown {
835552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
836552f7358SJed Brown   const PetscScalar x0        = vertices[0];
837552f7358SJed Brown   const PetscScalar y0        = vertices[1];
838552f7358SJed Brown   const PetscScalar z0        = vertices[2];
8397a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
8407a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
8417a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
842552f7358SJed Brown   const PetscScalar x2        = vertices[6];
843552f7358SJed Brown   const PetscScalar y2        = vertices[7];
844552f7358SJed Brown   const PetscScalar z2        = vertices[8];
8457a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
8467a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
8477a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
848552f7358SJed Brown   const PetscScalar x4        = vertices[12];
849552f7358SJed Brown   const PetscScalar y4        = vertices[13];
850552f7358SJed Brown   const PetscScalar z4        = vertices[14];
851552f7358SJed Brown   const PetscScalar x5        = vertices[15];
852552f7358SJed Brown   const PetscScalar y5        = vertices[16];
853552f7358SJed Brown   const PetscScalar z5        = vertices[17];
854552f7358SJed Brown   const PetscScalar x6        = vertices[18];
855552f7358SJed Brown   const PetscScalar y6        = vertices[19];
856552f7358SJed Brown   const PetscScalar z6        = vertices[20];
857552f7358SJed Brown   const PetscScalar x7        = vertices[21];
858552f7358SJed Brown   const PetscScalar y7        = vertices[22];
859552f7358SJed Brown   const PetscScalar z7        = vertices[23];
860552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
861552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
862552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
863552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
864552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
865552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
866552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
867552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
868552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
869552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
870552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
871552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
87256044e6dSMatthew G. Knepley   const PetscScalar *ref;
873552f7358SJed Brown   PetscErrorCode    ierr;
874552f7358SJed Brown 
875552f7358SJed Brown   PetscFunctionBegin;
87656044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
877552f7358SJed Brown   {
878552f7358SJed Brown     const PetscScalar x       = ref[0];
879552f7358SJed Brown     const PetscScalar y       = ref[1];
880552f7358SJed Brown     const PetscScalar z       = ref[2];
881552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
882da80777bSKarl Rupp     PetscScalar       values[9];
883da80777bSKarl Rupp 
884da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
885da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
886da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
887da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
888da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
889da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
890da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
891da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
892da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
8931aa26658SKarl Rupp 
89494ab13aaSBarry Smith     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
895552f7358SJed Brown   }
896552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
89756044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
89894ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89994ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
900552f7358SJed Brown   PetscFunctionReturn(0);
901552f7358SJed Brown }
902552f7358SJed Brown 
903a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
904a6dfd86eSKarl Rupp {
905fafc0619SMatthew G Knepley   DM             dmCoord;
906552f7358SJed Brown   SNES           snes;
907552f7358SJed Brown   KSP            ksp;
908552f7358SJed Brown   PC             pc;
909552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
910552f7358SJed Brown   Mat            J;
91156044e6dSMatthew G. Knepley   const PetscScalar *coords;
91256044e6dSMatthew G. Knepley   PetscScalar    *a;
913552f7358SJed Brown   PetscInt       p;
914552f7358SJed Brown   PetscErrorCode ierr;
915552f7358SJed Brown 
916552f7358SJed Brown   PetscFunctionBegin;
917552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
918fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
919552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
920552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
921552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
922552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
923c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
924552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
925552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
926552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
927552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
928552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
929552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
9300298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
9310298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
932552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
933552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
934552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
935552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
936552f7358SJed Brown 
93756044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
938552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
939552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
940a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
941552f7358SJed Brown     PetscScalar *xi;
942cb313848SJed Brown     PetscReal    xir[3];
943552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
944552f7358SJed Brown 
945552f7358SJed Brown     /* Can make this do all points at once */
9460298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
947ff1e0c32SBarry Smith     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
9480298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
949ff1e0c32SBarry Smith     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
9500298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
9510298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
952552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
953552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
954552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
955552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
956552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
957552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
958552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
959cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
960cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
961cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
962552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
963552f7358SJed Brown       a[p*ctx->dof+comp] =
964cb313848SJed Brown         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
9657a1931ceSMatthew G. Knepley         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
966cb313848SJed Brown         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
9677a1931ceSMatthew G. Knepley         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
968cb313848SJed Brown         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
969cb313848SJed Brown         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
970cb313848SJed Brown         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
971cb313848SJed Brown         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
972552f7358SJed Brown     }
973552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
9740298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
9750298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
976552f7358SJed Brown   }
977552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
97856044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
979552f7358SJed Brown 
980552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
981552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
982552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
983552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
984552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
985552f7358SJed Brown   PetscFunctionReturn(0);
986552f7358SJed Brown }
987552f7358SJed Brown 
9884267b1a3SMatthew G. Knepley /*@C
9894267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9904267b1a3SMatthew G. Knepley 
991552f7358SJed Brown   Input Parameters:
992552f7358SJed Brown + ctx - The DMInterpolationInfo context
993552f7358SJed Brown . dm  - The DM
994552f7358SJed Brown - x   - The local vector containing the field to be interpolated
995552f7358SJed Brown 
996552f7358SJed Brown   Output Parameters:
997552f7358SJed Brown . v   - The vector containing the interpolated values
9984267b1a3SMatthew G. Knepley 
9994267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
10004267b1a3SMatthew G. Knepley 
10014267b1a3SMatthew G. Knepley   Level: beginner
10024267b1a3SMatthew G. Knepley 
10034267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
10044267b1a3SMatthew G. Knepley @*/
10050adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
10060adebc6cSBarry Smith {
1007680d18d5SMatthew G. Knepley   PetscInt       n;
1008552f7358SJed Brown   PetscErrorCode ierr;
1009552f7358SJed Brown 
1010552f7358SJed Brown   PetscFunctionBegin;
1011552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1012552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1013552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
1014552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1015ff1e0c32SBarry Smith   if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %D should be %D", n, ctx->n*ctx->dof);
1016552f7358SJed Brown   if (n) {
1017680d18d5SMatthew G. Knepley     PetscDS        ds;
1018680d18d5SMatthew G. Knepley     DMPolytopeType ct;
1019680d18d5SMatthew G. Knepley     PetscBool      done = PETSC_FALSE;
1020680d18d5SMatthew G. Knepley 
1021680d18d5SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1022680d18d5SMatthew G. Knepley     if (ds) {
1023680d18d5SMatthew G. Knepley       const PetscScalar *coords;
1024680d18d5SMatthew G. Knepley       PetscScalar       *interpolant;
1025680d18d5SMatthew G. Knepley       PetscInt           cdim, d, p, Nf, field, c = 0;
1026680d18d5SMatthew G. Knepley 
1027680d18d5SMatthew G. Knepley       ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1028680d18d5SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
1029680d18d5SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
1030680d18d5SMatthew G. Knepley         PetscTabulation T;
1031680d18d5SMatthew G. Knepley         PetscFE         fe;
1032680d18d5SMatthew G. Knepley         PetscClassId    id;
1033680d18d5SMatthew G. Knepley         PetscReal       xi[3];
1034680d18d5SMatthew G. Knepley         PetscInt        Nc, f, fc;
1035680d18d5SMatthew G. Knepley 
1036680d18d5SMatthew G. Knepley         ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
1037680d18d5SMatthew G. Knepley         ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
1038680d18d5SMatthew G. Knepley         if (id != PETSCFE_CLASSID) break;
1039680d18d5SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1040680d18d5SMatthew G. Knepley         ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1041680d18d5SMatthew G. Knepley         ierr = VecGetArrayWrite(v, &interpolant);CHKERRQ(ierr);
1042680d18d5SMatthew G. Knepley         for (p = 0; p < ctx->n; ++p) {
1043680d18d5SMatthew G. Knepley           PetscScalar *xa = NULL;
1044680d18d5SMatthew G. Knepley           PetscReal    pcoords[3];
1045680d18d5SMatthew G. Knepley 
104652aa1562SMatthew G. Knepley           if (ctx->cells[p] < 0) continue;
1047680d18d5SMatthew G. Knepley           for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
1048680d18d5SMatthew G. Knepley           ierr = DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi);CHKERRQ(ierr);
1049680d18d5SMatthew G. Knepley           ierr = DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1050680d18d5SMatthew G. Knepley           ierr = PetscFECreateTabulation(fe, 1, 1, xi, 0, &T);CHKERRQ(ierr);
1051680d18d5SMatthew G. Knepley           {
1052680d18d5SMatthew G. Knepley             const PetscReal *basis = T->T[0];
1053680d18d5SMatthew G. Knepley             const PetscInt   Nb    = T->Nb;
1054680d18d5SMatthew G. Knepley             const PetscInt   Nc    = T->Nc;
1055680d18d5SMatthew G. Knepley             for (fc = 0; fc < Nc; ++fc) {
1056680d18d5SMatthew G. Knepley               interpolant[p*ctx->dof+c+fc] = 0.0;
1057680d18d5SMatthew G. Knepley               for (f = 0; f < Nb; ++f) {
1058680d18d5SMatthew G. Knepley                 interpolant[p*ctx->dof+c+fc] += xa[f]*basis[(0*Nb + f)*Nc + fc];
1059552f7358SJed Brown               }
1060680d18d5SMatthew G. Knepley             }
1061680d18d5SMatthew G. Knepley           }
1062680d18d5SMatthew G. Knepley           ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
1063680d18d5SMatthew G. Knepley           ierr = DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1064680d18d5SMatthew G. Knepley         }
1065680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayWrite(v, &interpolant);CHKERRQ(ierr);
1066680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1067680d18d5SMatthew G. Knepley         c += Nc;
1068680d18d5SMatthew G. Knepley       }
1069680d18d5SMatthew G. Knepley       if (field == Nf) {
1070680d18d5SMatthew G. Knepley         done = PETSC_TRUE;
1071680d18d5SMatthew G. Knepley         if (c != ctx->dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %D != %D dof specified for interpolation", c, ctx->dof);
1072680d18d5SMatthew G. Knepley       }
1073680d18d5SMatthew G. Knepley     }
1074680d18d5SMatthew G. Knepley     if (!done) {
1075680d18d5SMatthew G. Knepley       /* TODO Check each cell individually */
1076680d18d5SMatthew G. Knepley       ierr = DMPlexGetCellType(dm, ctx->cells[0], &ct);CHKERRQ(ierr);
1077680d18d5SMatthew G. Knepley       switch (ct) {
1078680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:      ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1079680d18d5SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL: ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1080680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TETRAHEDRON:   ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1081680d18d5SMatthew G. Knepley         case DM_POLYTOPE_HEXAHEDRON:    ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1082680d18d5SMatthew G. Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support fpr cell type %s", DMPolytopeTypes[ct]);
1083680d18d5SMatthew G. Knepley       }
1084680d18d5SMatthew G. Knepley     }
1085552f7358SJed Brown   }
1086552f7358SJed Brown   PetscFunctionReturn(0);
1087552f7358SJed Brown }
1088552f7358SJed Brown 
10894267b1a3SMatthew G. Knepley /*@C
10904267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
10914267b1a3SMatthew G. Knepley 
10924267b1a3SMatthew G. Knepley   Collective on ctx
10934267b1a3SMatthew G. Knepley 
10944267b1a3SMatthew G. Knepley   Input Parameter:
10954267b1a3SMatthew G. Knepley . ctx - the context
10964267b1a3SMatthew G. Knepley 
10974267b1a3SMatthew G. Knepley   Level: beginner
10984267b1a3SMatthew G. Knepley 
10994267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
11004267b1a3SMatthew G. Knepley @*/
11010adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
11020adebc6cSBarry Smith {
1103552f7358SJed Brown   PetscErrorCode ierr;
1104552f7358SJed Brown 
1105552f7358SJed Brown   PetscFunctionBegin;
1106552f7358SJed Brown   PetscValidPointer(ctx, 2);
1107552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
1108552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
1109552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
1110552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
11110298fd71SBarry Smith   *ctx = NULL;
1112552f7358SJed Brown   PetscFunctionReturn(0);
1113552f7358SJed Brown }
1114cc0c4584SMatthew G. Knepley 
1115cc0c4584SMatthew G. Knepley /*@C
1116cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1117cc0c4584SMatthew G. Knepley 
1118cc0c4584SMatthew G. Knepley   Collective on SNES
1119cc0c4584SMatthew G. Knepley 
1120cc0c4584SMatthew G. Knepley   Input Parameters:
1121cc0c4584SMatthew G. Knepley + snes   - the SNES context
1122cc0c4584SMatthew G. Knepley . its    - iteration number
1123cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1124d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
1125cc0c4584SMatthew G. Knepley 
1126cc0c4584SMatthew G. Knepley   Notes:
1127cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1128cc0c4584SMatthew G. Knepley 
1129cc0c4584SMatthew G. Knepley   Level: intermediate
1130cc0c4584SMatthew G. Knepley 
1131cc0c4584SMatthew G. Knepley .seealso: SNESMonitorSet(), SNESMonitorDefault()
1132cc0c4584SMatthew G. Knepley @*/
1133d43b4f6eSBarry Smith PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1134cc0c4584SMatthew G. Knepley {
1135d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1136cc0c4584SMatthew G. Knepley   Vec                res;
1137cc0c4584SMatthew G. Knepley   DM                 dm;
1138cc0c4584SMatthew G. Knepley   PetscSection       s;
1139cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1140cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1141cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1142cc0c4584SMatthew G. Knepley   PetscErrorCode     ierr;
1143cc0c4584SMatthew G. Knepley 
1144cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11454d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
11469e5d0892SLisandro Dalcin   ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr);
1147cc0c4584SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
114892fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
1149cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
1150cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1151cc0c4584SMatthew G. Knepley   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
1152cc0c4584SMatthew G. Knepley   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
1153cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1154cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1155cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1156cc0c4584SMatthew G. Knepley 
1157cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1158cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1159cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1160cc0c4584SMatthew G. Knepley     }
1161cc0c4584SMatthew G. Knepley   }
1162cc0c4584SMatthew G. Knepley   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
1163b2566f29SBarry Smith   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1164d43b4f6eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1165cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1166cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
1167cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
1168cc0c4584SMatthew G. Knepley     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
1169cc0c4584SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
1170cc0c4584SMatthew G. Knepley   }
1171cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
1172cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1173d43b4f6eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1174cc0c4584SMatthew G. Knepley   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
1175cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1176cc0c4584SMatthew G. Knepley }
117724cdb843SMatthew G. Knepley 
117824cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
117924cdb843SMatthew G. Knepley 
1180*6528b96dSMatthew G. Knepley PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
1181*6528b96dSMatthew G. Knepley {
1182*6528b96dSMatthew G. Knepley   PetscInt       depth;
1183*6528b96dSMatthew G. Knepley   PetscErrorCode ierr;
1184*6528b96dSMatthew G. Knepley 
1185*6528b96dSMatthew G. Knepley   PetscFunctionBegin;
1186*6528b96dSMatthew G. Knepley   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1187*6528b96dSMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, cellIS);CHKERRQ(ierr);
1188*6528b96dSMatthew G. Knepley   if (!*cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, cellIS);CHKERRQ(ierr);}
1189*6528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
1190*6528b96dSMatthew G. Knepley }
1191*6528b96dSMatthew G. Knepley 
119224cdb843SMatthew G. Knepley /*@
119324cdb843SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
119424cdb843SMatthew G. Knepley 
119524cdb843SMatthew G. Knepley   Input Parameters:
119624cdb843SMatthew G. Knepley + dm - The mesh
119724cdb843SMatthew G. Knepley . X  - Local solution
119824cdb843SMatthew G. Knepley - user - The user context
119924cdb843SMatthew G. Knepley 
120024cdb843SMatthew G. Knepley   Output Parameter:
120124cdb843SMatthew G. Knepley . F  - Local output vector
120224cdb843SMatthew G. Knepley 
120324cdb843SMatthew G. Knepley   Level: developer
120424cdb843SMatthew G. Knepley 
12057a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
120624cdb843SMatthew G. Knepley @*/
120724cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
120824cdb843SMatthew G. Knepley {
12096da023fcSToby Isaac   DM             plex;
1210083401c6SMatthew G. Knepley   IS             allcellIS;
1211*6528b96dSMatthew G. Knepley   PetscInt       Nds, s;
121224cdb843SMatthew G. Knepley   PetscErrorCode ierr;
121324cdb843SMatthew G. Knepley 
121424cdb843SMatthew G. Knepley   PetscFunctionBegin;
12156da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1216*6528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1217*6528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1218*6528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1219*6528b96dSMatthew G. Knepley     PetscDS          ds;
1220*6528b96dSMatthew G. Knepley     IS               cellIS;
1221*6528b96dSMatthew G. Knepley     PetscHashFormKey key;
1222*6528b96dSMatthew G. Knepley 
1223*6528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
1224*6528b96dSMatthew G. Knepley     key.value = 0;
1225*6528b96dSMatthew G. Knepley     key.field = 0;
1226*6528b96dSMatthew G. Knepley     if (!key.label) {
1227*6528b96dSMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1228*6528b96dSMatthew G. Knepley       cellIS = allcellIS;
1229*6528b96dSMatthew G. Knepley     } else {
1230*6528b96dSMatthew G. Knepley       IS pointIS;
1231*6528b96dSMatthew G. Knepley 
1232*6528b96dSMatthew G. Knepley       key.value = 1;
1233*6528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1234*6528b96dSMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1235*6528b96dSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1236*6528b96dSMatthew G. Knepley     }
1237*6528b96dSMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
1238*6528b96dSMatthew G. Knepley     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1239*6528b96dSMatthew G. Knepley   }
1240*6528b96dSMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
1241*6528b96dSMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1242*6528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
1243*6528b96dSMatthew G. Knepley }
1244*6528b96dSMatthew G. Knepley 
1245*6528b96dSMatthew G. Knepley PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user)
1246*6528b96dSMatthew G. Knepley {
1247*6528b96dSMatthew G. Knepley   DM             plex;
1248*6528b96dSMatthew G. Knepley   IS             allcellIS;
1249*6528b96dSMatthew G. Knepley   PetscInt       Nds, s;
1250*6528b96dSMatthew G. Knepley   PetscErrorCode ierr;
1251*6528b96dSMatthew G. Knepley 
1252*6528b96dSMatthew G. Knepley   PetscFunctionBegin;
1253*6528b96dSMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1254*6528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1255*6528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1256083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1257083401c6SMatthew G. Knepley     PetscDS ds;
1258083401c6SMatthew G. Knepley     DMLabel label;
1259083401c6SMatthew G. Knepley     IS      cellIS;
1260083401c6SMatthew G. Knepley 
1261083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1262*6528b96dSMatthew G. Knepley     {
1263*6528b96dSMatthew G. Knepley       PetscHMapForm     resmap[2] = {ds->wf->f0, ds->wf->f1};
1264*6528b96dSMatthew G. Knepley       PetscWeakForm     wf;
1265*6528b96dSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
1266*6528b96dSMatthew G. Knepley       PetscHashFormKey *reskeys;
1267*6528b96dSMatthew G. Knepley 
1268*6528b96dSMatthew G. Knepley       /* Get unique residual keys */
1269*6528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
1270*6528b96dSMatthew G. Knepley         PetscInt Nkm;
1271*6528b96dSMatthew G. Knepley         ierr = PetscHMapFormGetSize(resmap[m], &Nkm);CHKERRQ(ierr);
1272*6528b96dSMatthew G. Knepley         Nk  += Nkm;
1273*6528b96dSMatthew G. Knepley       }
1274*6528b96dSMatthew G. Knepley       ierr = PetscMalloc1(Nk, &reskeys);CHKERRQ(ierr);
1275*6528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
1276*6528b96dSMatthew G. Knepley         ierr = PetscHMapFormGetKeys(resmap[m], &off, reskeys);CHKERRQ(ierr);
1277*6528b96dSMatthew G. Knepley       }
1278*6528b96dSMatthew G. Knepley       if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1279*6528b96dSMatthew G. Knepley       ierr = PetscHashFormKeySort(Nk, reskeys);CHKERRQ(ierr);
1280*6528b96dSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
1281*6528b96dSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
1282*6528b96dSMatthew G. Knepley           ++k;
1283*6528b96dSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
1284*6528b96dSMatthew G. Knepley         }
1285*6528b96dSMatthew G. Knepley       }
1286*6528b96dSMatthew G. Knepley       Nk = k;
1287*6528b96dSMatthew G. Knepley 
1288*6528b96dSMatthew G. Knepley       ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1289*6528b96dSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
1290*6528b96dSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
1291*6528b96dSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
1292*6528b96dSMatthew G. Knepley 
1293083401c6SMatthew G. Knepley         if (!label) {
1294083401c6SMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1295083401c6SMatthew G. Knepley           cellIS = allcellIS;
1296083401c6SMatthew G. Knepley         } else {
1297083401c6SMatthew G. Knepley           IS pointIS;
1298083401c6SMatthew G. Knepley 
1299*6528b96dSMatthew G. Knepley           ierr = DMLabelGetStratumIS(label, val, &pointIS);CHKERRQ(ierr);
1300083401c6SMatthew G. Knepley           ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1301083401c6SMatthew G. Knepley           ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
13024a3e9fdbSToby Isaac         }
1303*6528b96dSMatthew G. Knepley         ierr = DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
13044a3e9fdbSToby Isaac         ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1305083401c6SMatthew G. Knepley       }
1306*6528b96dSMatthew G. Knepley       ierr = PetscFree(reskeys);CHKERRQ(ierr);
1307*6528b96dSMatthew G. Knepley     }
1308*6528b96dSMatthew G. Knepley   }
1309083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
13109a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
131124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
131224cdb843SMatthew G. Knepley }
131324cdb843SMatthew G. Knepley 
1314bdd6f66aSToby Isaac /*@
1315bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1316bdd6f66aSToby Isaac 
1317bdd6f66aSToby Isaac   Input Parameters:
1318bdd6f66aSToby Isaac + dm - The mesh
1319bdd6f66aSToby Isaac - user - The user context
1320bdd6f66aSToby Isaac 
1321bdd6f66aSToby Isaac   Output Parameter:
1322bdd6f66aSToby Isaac . X  - Local solution
1323bdd6f66aSToby Isaac 
1324bdd6f66aSToby Isaac   Level: developer
1325bdd6f66aSToby Isaac 
13267a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
1327bdd6f66aSToby Isaac @*/
1328bdd6f66aSToby Isaac PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1329bdd6f66aSToby Isaac {
1330bdd6f66aSToby Isaac   DM             plex;
1331bdd6f66aSToby Isaac   PetscErrorCode ierr;
1332bdd6f66aSToby Isaac 
1333bdd6f66aSToby Isaac   PetscFunctionBegin;
1334bdd6f66aSToby Isaac   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1335bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1336bdd6f66aSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1337bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1338bdd6f66aSToby Isaac }
1339bdd6f66aSToby Isaac 
13407a73cf09SMatthew G. Knepley /*@
13417a73cf09SMatthew G. Knepley   DMPlexComputeJacobianAction - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.
13427a73cf09SMatthew G. Knepley 
13437a73cf09SMatthew G. Knepley   Input Parameters:
13447a73cf09SMatthew G. Knepley + dm - The mesh
13457a73cf09SMatthew G. Knepley . cellIS -
13467a73cf09SMatthew G. Knepley . t  - The time
13477a73cf09SMatthew G. Knepley . X_tShift - The multiplier for the Jacobian with repsect to X_t
13487a73cf09SMatthew G. Knepley . X  - Local solution vector
13497a73cf09SMatthew G. Knepley . X_t  - Time-derivative of the local solution vector
13507a73cf09SMatthew G. Knepley . Y  - Local input vector
13517a73cf09SMatthew G. Knepley - user - The user context
13527a73cf09SMatthew G. Knepley 
13537a73cf09SMatthew G. Knepley   Output Parameter:
13547a73cf09SMatthew G. Knepley . Z - Local output vector
13557a73cf09SMatthew G. Knepley 
13567a73cf09SMatthew G. Knepley   Note:
13577a73cf09SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
13587a73cf09SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
13597a73cf09SMatthew G. Knepley 
13607a73cf09SMatthew G. Knepley   Level: developer
13617a73cf09SMatthew G. Knepley 
13627a73cf09SMatthew G. Knepley .seealso: FormFunctionLocal()
13637a73cf09SMatthew G. Knepley @*/
13647a73cf09SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
1365a925c78cSMatthew G. Knepley {
1366a925c78cSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1367a925c78cSMatthew G. Knepley   const char       *name  = "Jacobian";
13687a73cf09SMatthew G. Knepley   DM                dmAux, plex, plexAux = NULL;
1369a6e0b375SMatthew G. Knepley   DMEnclosureType   encAux;
1370c330f8ffSToby Isaac   Vec               A;
1371a925c78cSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
1372a925c78cSMatthew G. Knepley   PetscQuadrature   quad;
1373a925c78cSMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1374a925c78cSMatthew G. Knepley   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
13759044fa66SMatthew G. Knepley   PetscInt          Nf, fieldI, fieldJ;
13764d0b9603SSander Arens   PetscInt          totDim, totDimAux = 0;
13779044fa66SMatthew G. Knepley   const PetscInt   *cells;
13789044fa66SMatthew G. Knepley   PetscInt          cStart, cEnd, numCells, c;
137975b37a90SMatthew G. Knepley   PetscBool         hasDyn;
1380c330f8ffSToby Isaac   DMField           coordField;
1381*6528b96dSMatthew G. Knepley   PetscHashFormKey  key;
1382a925c78cSMatthew G. Knepley   PetscErrorCode    ierr;
1383a925c78cSMatthew G. Knepley 
1384a925c78cSMatthew G. Knepley   PetscFunctionBegin;
1385a925c78cSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
13867a73cf09SMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
13877a73cf09SMatthew G. Knepley   if (!cellIS) {
13887a73cf09SMatthew G. Knepley     PetscInt depth;
13897a73cf09SMatthew G. Knepley 
13907a73cf09SMatthew G. Knepley     ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
13917a73cf09SMatthew G. Knepley     ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
13927a73cf09SMatthew G. Knepley     if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
13937a73cf09SMatthew G. Knepley   } else {
13947a73cf09SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) cellIS);CHKERRQ(ierr);
13957a73cf09SMatthew G. Knepley   }
1396*6528b96dSMatthew G. Knepley   key.label = NULL;
1397*6528b96dSMatthew G. Knepley   key.value = 0;
139892fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
1399e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1400d28747ceSMatthew G. Knepley   ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr);
1401d28747ceSMatthew G. Knepley   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1402d28747ceSMatthew G. Knepley   ierr = DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);CHKERRQ(ierr);
1403d28747ceSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1404a925c78cSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1405a925c78cSMatthew G. Knepley   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
1406a925c78cSMatthew G. Knepley   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1407a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1408a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1409a925c78cSMatthew G. Knepley   if (dmAux) {
1410a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
14117a73cf09SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
141292fd8e1eSJed Brown     ierr = DMGetLocalSection(plexAux, &sectionAux);CHKERRQ(ierr);
1413a925c78cSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1414a925c78cSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1415a925c78cSMatthew G. Knepley   }
1416a925c78cSMatthew G. Knepley   ierr = VecSet(Z, 0.0);CHKERRQ(ierr);
1417a925c78cSMatthew G. Knepley   ierr = PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);CHKERRQ(ierr);
1418a925c78cSMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
14194a3e9fdbSToby Isaac   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1420a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
14219044fa66SMatthew G. Knepley     const PetscInt cell = cells ? cells[c] : c;
14229044fa66SMatthew G. Knepley     const PetscInt cind = c - cStart;
1423a925c78cSMatthew G. Knepley     PetscScalar   *x = NULL,  *x_t = NULL;
1424a925c78cSMatthew G. Knepley     PetscInt       i;
1425a925c78cSMatthew G. Knepley 
14269044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
14279044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
14289044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
1429a925c78cSMatthew G. Knepley     if (X_t) {
14309044fa66SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
14319044fa66SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
14329044fa66SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
1433a925c78cSMatthew G. Knepley     }
1434a925c78cSMatthew G. Knepley     if (dmAux) {
143544171101SMatthew G. Knepley       PetscInt subcell;
1436a6e0b375SMatthew G. Knepley       ierr = DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);CHKERRQ(ierr);
143744171101SMatthew G. Knepley       ierr = DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
14389044fa66SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
143944171101SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
1440a925c78cSMatthew G. Knepley     }
14419044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
14429044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
14439044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
1444a925c78cSMatthew G. Knepley   }
1445580bdb30SBarry Smith   ierr = PetscArrayzero(elemMat, numCells*totDim*totDim);CHKERRQ(ierr);
1446580bdb30SBarry Smith   if (hasDyn)  {ierr = PetscArrayzero(elemMatD, numCells*totDim*totDim);CHKERRQ(ierr);}
1447a925c78cSMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1448a925c78cSMatthew G. Knepley     PetscFE  fe;
1449c330f8ffSToby Isaac     PetscInt Nb;
1450a925c78cSMatthew G. Knepley     /* Conforming batches */
1451a925c78cSMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1452a925c78cSMatthew G. Knepley     /* Remainder */
1453c330f8ffSToby Isaac     PetscInt Nr, offset, Nq;
1454c330f8ffSToby Isaac     PetscQuadrature qGeom = NULL;
1455b7260050SToby Isaac     PetscInt    maxDegree;
1456c330f8ffSToby Isaac     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
1457a925c78cSMatthew G. Knepley 
1458a925c78cSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1459a925c78cSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1460a925c78cSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1461a925c78cSMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1462b7260050SToby Isaac     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
1463b7260050SToby Isaac     if (maxDegree <= 1) {ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);}
1464c330f8ffSToby Isaac     if (!qGeom) {
1465c330f8ffSToby Isaac       ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr);
1466c330f8ffSToby Isaac       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1467c330f8ffSToby Isaac     }
1468c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
14694a3e9fdbSToby Isaac     ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1470c330f8ffSToby Isaac     blockSize = Nb;
1471a925c78cSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1472a925c78cSMatthew G. Knepley     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1473a925c78cSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1474a925c78cSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1475a925c78cSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1476a925c78cSMatthew G. Knepley     offset    = numCells - Nr;
1477c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
1478c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1479a925c78cSMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1480*6528b96dSMatthew G. Knepley       key.field = fieldI*Nf + fieldJ;
1481*6528b96dSMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
1482*6528b96dSMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
1483a925c78cSMatthew G. Knepley       if (hasDyn) {
1484*6528b96dSMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
1485*6528b96dSMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr);
1486a925c78cSMatthew G. Knepley       }
1487a925c78cSMatthew G. Knepley     }
1488c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1489c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
14904a3e9fdbSToby Isaac     ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1491c330f8ffSToby Isaac     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
1492a925c78cSMatthew G. Knepley   }
1493a925c78cSMatthew G. Knepley   if (hasDyn) {
14949044fa66SMatthew G. Knepley     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
1495a925c78cSMatthew G. Knepley   }
1496a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
14979044fa66SMatthew G. Knepley     const PetscInt     cell = cells ? cells[c] : c;
14989044fa66SMatthew G. Knepley     const PetscInt     cind = c - cStart;
1499a925c78cSMatthew G. Knepley     const PetscBLASInt M = totDim, one = 1;
1500a925c78cSMatthew G. Knepley     const PetscScalar  a = 1.0, b = 0.0;
1501a925c78cSMatthew G. Knepley 
15029044fa66SMatthew G. Knepley     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
1503a925c78cSMatthew G. Knepley     if (mesh->printFEM > 1) {
15049044fa66SMatthew G. Knepley       ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);
15059044fa66SMatthew G. Knepley       ierr = DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);CHKERRQ(ierr);
1506a925c78cSMatthew G. Knepley       ierr = DMPrintCellVector(c, "Z",  totDim, z);CHKERRQ(ierr);
1507a925c78cSMatthew G. Knepley     }
15089044fa66SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);CHKERRQ(ierr);
1509a925c78cSMatthew G. Knepley   }
1510a925c78cSMatthew G. Knepley   ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr);
1511a925c78cSMatthew G. Knepley   if (mesh->printFEM) {
1512ea13f565SStefano Zampini     ierr = PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n");CHKERRQ(ierr);
1513a8c5bd36SStefano Zampini     ierr = VecView(Z, NULL);CHKERRQ(ierr);
1514a925c78cSMatthew G. Knepley   }
15157a73cf09SMatthew G. Knepley   ierr = PetscFree(a);CHKERRQ(ierr);
15167a73cf09SMatthew G. Knepley   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
15177a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plexAux);CHKERRQ(ierr);
15187a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1519a925c78cSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1520a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1521a925c78cSMatthew G. Knepley }
1522a925c78cSMatthew G. Knepley 
152324cdb843SMatthew G. Knepley /*@
152424cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
152524cdb843SMatthew G. Knepley 
152624cdb843SMatthew G. Knepley   Input Parameters:
152724cdb843SMatthew G. Knepley + dm - The mesh
152824cdb843SMatthew G. Knepley . X  - Local input vector
152924cdb843SMatthew G. Knepley - user - The user context
153024cdb843SMatthew G. Knepley 
153124cdb843SMatthew G. Knepley   Output Parameter:
153224cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
153324cdb843SMatthew G. Knepley 
153424cdb843SMatthew G. Knepley   Note:
153524cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
153624cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
153724cdb843SMatthew G. Knepley 
153824cdb843SMatthew G. Knepley   Level: developer
153924cdb843SMatthew G. Knepley 
154024cdb843SMatthew G. Knepley .seealso: FormFunctionLocal()
154124cdb843SMatthew G. Knepley @*/
154224cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
154324cdb843SMatthew G. Knepley {
15446da023fcSToby Isaac   DM             plex;
1545083401c6SMatthew G. Knepley   IS             allcellIS;
1546f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
1547*6528b96dSMatthew G. Knepley   PetscInt       Nds, s;
154824cdb843SMatthew G. Knepley   PetscErrorCode ierr;
154924cdb843SMatthew G. Knepley 
155024cdb843SMatthew G. Knepley   PetscFunctionBegin;
15516da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1552*6528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1553*6528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1554083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1555083401c6SMatthew G. Knepley     PetscDS          ds;
1556083401c6SMatthew G. Knepley     IS               cellIS;
1557*6528b96dSMatthew G. Knepley     PetscHashFormKey key;
1558083401c6SMatthew G. Knepley 
1559*6528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
1560*6528b96dSMatthew G. Knepley     key.value = 0;
1561*6528b96dSMatthew G. Knepley     key.field = 0;
1562*6528b96dSMatthew G. Knepley     if (!key.label) {
1563083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1564083401c6SMatthew G. Knepley       cellIS = allcellIS;
1565083401c6SMatthew G. Knepley     } else {
1566083401c6SMatthew G. Knepley       IS pointIS;
1567083401c6SMatthew G. Knepley 
1568*6528b96dSMatthew G. Knepley       key.value = 1;
1569*6528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1570083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1571083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1572083401c6SMatthew G. Knepley     }
1573083401c6SMatthew G. Knepley     if (!s) {
1574083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1575083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1576f04eb4edSMatthew G. Knepley       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
1577f04eb4edSMatthew G. Knepley       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1578083401c6SMatthew G. Knepley     }
1579*6528b96dSMatthew G. Knepley     ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
15804a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1581083401c6SMatthew G. Knepley   }
1582083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
15839a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
158424cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
158524cdb843SMatthew G. Knepley }
15861878804aSMatthew G. Knepley 
158738cfc46eSPierre Jolivet /*
158838cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
158938cfc46eSPierre Jolivet 
159038cfc46eSPierre Jolivet    Input Parameters:
159138cfc46eSPierre Jolivet +     X - SNES linearization point
159238cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
159338cfc46eSPierre Jolivet 
159438cfc46eSPierre Jolivet    Output Parameter:
159538cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
159638cfc46eSPierre Jolivet 
159738cfc46eSPierre Jolivet    Level: intermediate
159838cfc46eSPierre Jolivet 
159938cfc46eSPierre Jolivet .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
160038cfc46eSPierre Jolivet */
160138cfc46eSPierre Jolivet static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
160238cfc46eSPierre Jolivet {
160338cfc46eSPierre Jolivet   SNES           snes;
160438cfc46eSPierre Jolivet   Mat            pJ;
160538cfc46eSPierre Jolivet   DM             ovldm,origdm;
160638cfc46eSPierre Jolivet   DMSNES         sdm;
160738cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM,Vec,void*);
160838cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
160938cfc46eSPierre Jolivet   void           *bctx,*jctx;
161038cfc46eSPierre Jolivet   PetscErrorCode ierr;
161138cfc46eSPierre Jolivet 
161238cfc46eSPierre Jolivet   PetscFunctionBegin;
161338cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr);
161438cfc46eSPierre Jolivet   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr);
161538cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr);
161638cfc46eSPierre Jolivet   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr);
161738cfc46eSPierre Jolivet   ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr);
161838cfc46eSPierre Jolivet   ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr);
161938cfc46eSPierre Jolivet   ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr);
162038cfc46eSPierre Jolivet   ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr);
162138cfc46eSPierre Jolivet   ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr);
162238cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr);
162338cfc46eSPierre Jolivet   if (!snes) {
162438cfc46eSPierre Jolivet     ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr);
162538cfc46eSPierre Jolivet     ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr);
162638cfc46eSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr);
162738cfc46eSPierre Jolivet     ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr);
162838cfc46eSPierre Jolivet   }
162938cfc46eSPierre Jolivet   ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr);
163038cfc46eSPierre Jolivet   ierr = VecLockReadPush(X);CHKERRQ(ierr);
163138cfc46eSPierre Jolivet   PetscStackPush("SNES user Jacobian function");
163238cfc46eSPierre Jolivet   ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr);
163338cfc46eSPierre Jolivet   PetscStackPop;
163438cfc46eSPierre Jolivet   ierr = VecLockReadPop(X);CHKERRQ(ierr);
163538cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
163638cfc46eSPierre Jolivet   {
163738cfc46eSPierre Jolivet     Mat locpJ;
163838cfc46eSPierre Jolivet 
163938cfc46eSPierre Jolivet     ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr);
164038cfc46eSPierre Jolivet     ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
164138cfc46eSPierre Jolivet   }
164238cfc46eSPierre Jolivet   PetscFunctionReturn(0);
164338cfc46eSPierre Jolivet }
164438cfc46eSPierre Jolivet 
1645a925c78cSMatthew G. Knepley /*@
16469f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
16479f520fc2SToby Isaac 
16489f520fc2SToby Isaac   Input Parameters:
16499f520fc2SToby Isaac + dm - The DM object
1650dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
16519f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
16529f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
16531a244344SSatish Balay 
16541a244344SSatish Balay   Level: developer
16559f520fc2SToby Isaac @*/
16569f520fc2SToby Isaac PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
16579f520fc2SToby Isaac {
16589f520fc2SToby Isaac   PetscErrorCode ierr;
16599f520fc2SToby Isaac 
16609f520fc2SToby Isaac   PetscFunctionBegin;
16619f520fc2SToby Isaac   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
16629f520fc2SToby Isaac   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
16639f520fc2SToby Isaac   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
166438cfc46eSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr);
16659f520fc2SToby Isaac   PetscFunctionReturn(0);
16669f520fc2SToby Isaac }
16679f520fc2SToby Isaac 
1668f017bcb6SMatthew G. Knepley /*@C
1669f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1670f017bcb6SMatthew G. Knepley 
1671f017bcb6SMatthew G. Knepley   Input Parameters:
1672f017bcb6SMatthew G. Knepley + snes - the SNES object
1673f017bcb6SMatthew G. Knepley . dm   - the DM
1674f2cacb80SMatthew G. Knepley . t    - the time
1675f017bcb6SMatthew G. Knepley . u    - a DM vector
1676f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1677f017bcb6SMatthew G. Knepley 
1678f3c94560SMatthew G. Knepley   Output Parameters:
1679f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1680f3c94560SMatthew G. Knepley 
16817f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
16827f96f943SMatthew G. Knepley 
1683f017bcb6SMatthew G. Knepley   Level: developer
1684f017bcb6SMatthew G. Knepley 
16857f96f943SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1686f017bcb6SMatthew G. Knepley @*/
1687f2cacb80SMatthew G. Knepley PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
16881878804aSMatthew G. Knepley {
1689f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1690f017bcb6SMatthew G. Knepley   void            **ectxs;
1691f3c94560SMatthew G. Knepley   PetscReal        *err;
16927f96f943SMatthew G. Knepley   MPI_Comm          comm;
16937f96f943SMatthew G. Knepley   PetscInt          Nf, f;
16941878804aSMatthew G. Knepley   PetscErrorCode    ierr;
16951878804aSMatthew G. Knepley 
16961878804aSMatthew G. Knepley   PetscFunctionBegin;
1697f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1698f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1699f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1700534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
17017f96f943SMatthew G. Knepley 
1702f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
17037f96f943SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
17047f96f943SMatthew G. Knepley 
1705f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1706f017bcb6SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1707083401c6SMatthew G. Knepley   ierr = PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr);
17087f96f943SMatthew G. Knepley   {
17097f96f943SMatthew G. Knepley     PetscInt Nds, s;
17107f96f943SMatthew G. Knepley 
1711083401c6SMatthew G. Knepley     ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1712083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
17137f96f943SMatthew G. Knepley       PetscDS         ds;
1714083401c6SMatthew G. Knepley       DMLabel         label;
1715083401c6SMatthew G. Knepley       IS              fieldIS;
17167f96f943SMatthew G. Knepley       const PetscInt *fields;
17177f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1718083401c6SMatthew G. Knepley 
1719083401c6SMatthew G. Knepley       ierr = DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
1720083401c6SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
1721083401c6SMatthew G. Knepley       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
1722083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1723083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
1724083401c6SMatthew G. Knepley         ierr = PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);CHKERRQ(ierr);
1725083401c6SMatthew G. Knepley       }
1726083401c6SMatthew G. Knepley       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
1727083401c6SMatthew G. Knepley     }
1728083401c6SMatthew G. Knepley   }
1729f017bcb6SMatthew G. Knepley   if (Nf > 1) {
1730f2cacb80SMatthew G. Knepley     ierr = DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);CHKERRQ(ierr);
1731f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1732f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1733f3c94560SMatthew G. Knepley         if (err[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) err[f], f, (double) tol);
17341878804aSMatthew G. Knepley       }
1735b878532fSMatthew G. Knepley     } else if (error) {
1736b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
17371878804aSMatthew G. Knepley     } else {
1738f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
1739f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1740f017bcb6SMatthew G. Knepley         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
1741f3c94560SMatthew G. Knepley         ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr);
17421878804aSMatthew G. Knepley       }
1743f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
1744f017bcb6SMatthew G. Knepley     }
1745f017bcb6SMatthew G. Knepley   } else {
1746f2cacb80SMatthew G. Knepley     ierr = DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);CHKERRQ(ierr);
1747f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1748f3c94560SMatthew G. Knepley       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1749b878532fSMatthew G. Knepley     } else if (error) {
1750b878532fSMatthew G. Knepley       error[0] = err[0];
1751f017bcb6SMatthew G. Knepley     } else {
1752f3c94560SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);CHKERRQ(ierr);
1753f017bcb6SMatthew G. Knepley     }
1754f017bcb6SMatthew G. Knepley   }
1755f3c94560SMatthew G. Knepley   ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr);
1756f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1757f017bcb6SMatthew G. Knepley }
1758f017bcb6SMatthew G. Knepley 
1759f017bcb6SMatthew G. Knepley /*@C
1760f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1761f017bcb6SMatthew G. Knepley 
1762f017bcb6SMatthew G. Knepley   Input Parameters:
1763f017bcb6SMatthew G. Knepley + snes - the SNES object
1764f017bcb6SMatthew G. Knepley . dm   - the DM
1765f017bcb6SMatthew G. Knepley . u    - a DM vector
1766f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1767f017bcb6SMatthew G. Knepley 
1768f3c94560SMatthew G. Knepley   Output Parameters:
1769f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1770f3c94560SMatthew G. Knepley 
1771f017bcb6SMatthew G. Knepley   Level: developer
1772f017bcb6SMatthew G. Knepley 
1773f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1774f017bcb6SMatthew G. Knepley @*/
1775f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1776f017bcb6SMatthew G. Knepley {
1777f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1778f017bcb6SMatthew G. Knepley   Vec            r;
1779f017bcb6SMatthew G. Knepley   PetscReal      res;
1780f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1781f017bcb6SMatthew G. Knepley 
1782f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1783f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1784f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1785f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1786534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
1787f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1788f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1789f017bcb6SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
17901878804aSMatthew G. Knepley   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
17911878804aSMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
1792f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
1793f017bcb6SMatthew G. Knepley     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1794b878532fSMatthew G. Knepley   } else if (residual) {
1795b878532fSMatthew G. Knepley     *residual = res;
1796f017bcb6SMatthew G. Knepley   } else {
1797f017bcb6SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
17981878804aSMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
17991878804aSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
1800685405a1SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
1801685405a1SBarry Smith     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
1802f017bcb6SMatthew G. Knepley   }
1803f017bcb6SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
1804f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1805f017bcb6SMatthew G. Knepley }
1806f017bcb6SMatthew G. Knepley 
1807f017bcb6SMatthew G. Knepley /*@C
1808f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1809f017bcb6SMatthew G. Knepley 
1810f017bcb6SMatthew G. Knepley   Input Parameters:
1811f017bcb6SMatthew G. Knepley + snes - the SNES object
1812f017bcb6SMatthew G. Knepley . dm   - the DM
1813f017bcb6SMatthew G. Knepley . u    - a DM vector
1814f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1815f017bcb6SMatthew G. Knepley 
1816f3c94560SMatthew G. Knepley   Output Parameters:
1817f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1818f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1819f3c94560SMatthew G. Knepley 
1820f017bcb6SMatthew G. Knepley   Level: developer
1821f017bcb6SMatthew G. Knepley 
1822f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1823f017bcb6SMatthew G. Knepley @*/
1824f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1825f017bcb6SMatthew G. Knepley {
1826f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1827f017bcb6SMatthew G. Knepley   PetscDS        ds;
1828f017bcb6SMatthew G. Knepley   Mat            J, M;
1829f017bcb6SMatthew G. Knepley   MatNullSpace   nullspace;
1830f3c94560SMatthew G. Knepley   PetscReal      slope, intercept;
18317f96f943SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1832f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1833f017bcb6SMatthew G. Knepley 
1834f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1835f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1836f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1837f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1838534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1839534a8f05SLisandro Dalcin   if (convRate) PetscValidRealPointer(convRate, 5);
1840f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1841f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1842f017bcb6SMatthew G. Knepley   /* Create and view matrices */
1843f017bcb6SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
18447f96f943SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1845f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1846f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1847282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
1848282e7bb4SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
1849282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
1850f017bcb6SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
1851282e7bb4SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
1852282e7bb4SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
1853282e7bb4SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
1854282e7bb4SMatthew G. Knepley   } else {
1855282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
1856282e7bb4SMatthew G. Knepley   }
1857f017bcb6SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
1858282e7bb4SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
1859282e7bb4SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
1860f017bcb6SMatthew G. Knepley   /* Check nullspace */
1861f017bcb6SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
1862f017bcb6SMatthew G. Knepley   if (nullspace) {
18631878804aSMatthew G. Knepley     PetscBool isNull;
1864f017bcb6SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
1865f017bcb6SMatthew G. Knepley     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
18661878804aSMatthew G. Knepley   }
1867f017bcb6SMatthew G. Knepley   /* Taylor test */
1868f017bcb6SMatthew G. Knepley   {
1869f017bcb6SMatthew G. Knepley     PetscRandom rand;
1870f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1871f3c94560SMatthew G. Knepley     PetscReal   h;
1872f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1873f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1874f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1875f017bcb6SMatthew G. Knepley 
1876f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
1877f017bcb6SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
1878f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
1879f017bcb6SMatthew G. Knepley     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
1880f017bcb6SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
1881f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
1882f017bcb6SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
1883f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
1884f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1885f017bcb6SMatthew G. Knepley     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1886f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1887f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1888f017bcb6SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
1889f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
1890f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
1891f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1892f017bcb6SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
1893f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1894f017bcb6SMatthew G. Knepley       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
1895f017bcb6SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
1896f017bcb6SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
1897f017bcb6SMatthew G. Knepley 
1898f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1899f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1900f017bcb6SMatthew G. Knepley     }
1901f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
1902f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
1903f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
19041878804aSMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
1905f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
1906f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1907f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1908f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
1909f017bcb6SMatthew G. Knepley     }
1910f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
1911f017bcb6SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
1912f017bcb6SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
1913f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1914f017bcb6SMatthew G. Knepley     if (tol >= 0) {
1915b878532fSMatthew G. Knepley       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1916b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1917b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1918b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1919f017bcb6SMatthew G. Knepley     } else {
1920f3c94560SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
1921f017bcb6SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
1922f017bcb6SMatthew G. Knepley     }
1923f017bcb6SMatthew G. Knepley   }
19241878804aSMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
19251878804aSMatthew G. Knepley   PetscFunctionReturn(0);
19261878804aSMatthew G. Knepley }
19271878804aSMatthew G. Knepley 
19287f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1929f017bcb6SMatthew G. Knepley {
1930f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1931f017bcb6SMatthew G. Knepley 
1932f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1933f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);CHKERRQ(ierr);
1934f3c94560SMatthew G. Knepley   ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr);
1935f3c94560SMatthew G. Knepley   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr);
1936f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1937f017bcb6SMatthew G. Knepley }
1938f017bcb6SMatthew G. Knepley 
1939bee9a294SMatthew G. Knepley /*@C
1940bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1941bee9a294SMatthew G. Knepley 
1942bee9a294SMatthew G. Knepley   Input Parameters:
1943bee9a294SMatthew G. Knepley + snes - the SNES object
19447f96f943SMatthew G. Knepley - u    - representative SNES vector
19457f96f943SMatthew G. Knepley 
19467f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1947bee9a294SMatthew G. Knepley 
1948bee9a294SMatthew G. Knepley   Level: developer
1949bee9a294SMatthew G. Knepley @*/
19507f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
19511878804aSMatthew G. Knepley {
19521878804aSMatthew G. Knepley   DM             dm;
19531878804aSMatthew G. Knepley   Vec            sol;
19541878804aSMatthew G. Knepley   PetscBool      check;
19551878804aSMatthew G. Knepley   PetscErrorCode ierr;
19561878804aSMatthew G. Knepley 
19571878804aSMatthew G. Knepley   PetscFunctionBegin;
1958c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
19591878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
19601878804aSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
19611878804aSMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
19621878804aSMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
19637f96f943SMatthew G. Knepley   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
19641878804aSMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1965552f7358SJed Brown   PetscFunctionReturn(0);
1966552f7358SJed Brown }
1967