xref: /petsc/src/snes/utils/dmplexsnes.c (revision 8e3a2eeff93b8092aa720bf50e8b48bf277d7f13)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>     /*I "petscsnes.h"   I*/
324cdb843SMatthew G. Knepley #include <petscds.h>
4af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
5af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
6552f7358SJed Brown 
7649ef022SMatthew Knepley static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
8649ef022SMatthew Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9649ef022SMatthew Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10649ef022SMatthew Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
11649ef022SMatthew Knepley {
12649ef022SMatthew Knepley   p[0] = u[uOff[1]];
13649ef022SMatthew Knepley }
14649ef022SMatthew Knepley 
15649ef022SMatthew Knepley /*
16649ef022SMatthew 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.
17649ef022SMatthew Knepley 
18649ef022SMatthew Knepley   Collective on SNES
19649ef022SMatthew Knepley 
20649ef022SMatthew Knepley   Input Parameters:
21649ef022SMatthew Knepley + snes      - The SNES
22649ef022SMatthew Knepley . pfield    - The field number for pressure
23649ef022SMatthew Knepley . nullspace - The pressure nullspace
24649ef022SMatthew Knepley . u         - The solution vector
25649ef022SMatthew Knepley - ctx       - An optional user context
26649ef022SMatthew Knepley 
27649ef022SMatthew Knepley   Output Parameter:
28649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
29649ef022SMatthew Knepley 
30649ef022SMatthew Knepley   Notes:
31649ef022SMatthew 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.
32649ef022SMatthew Knepley 
33649ef022SMatthew Knepley   Level: developer
34649ef022SMatthew Knepley 
35649ef022SMatthew Knepley .seealso: SNESConvergedCorrectPressure()
36649ef022SMatthew Knepley */
37649ef022SMatthew Knepley static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
38649ef022SMatthew Knepley {
39649ef022SMatthew Knepley   DM             dm;
40649ef022SMatthew Knepley   PetscDS        ds;
41649ef022SMatthew Knepley   const Vec     *nullvecs;
42649ef022SMatthew Knepley   PetscScalar    pintd, *intc, *intn;
43649ef022SMatthew Knepley   MPI_Comm       comm;
44649ef022SMatthew Knepley   PetscInt       Nf, Nv;
45649ef022SMatthew Knepley   PetscErrorCode ierr;
46649ef022SMatthew Knepley 
47649ef022SMatthew Knepley   PetscFunctionBegin;
48649ef022SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
49649ef022SMatthew Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
50649ef022SMatthew Knepley   if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
51649ef022SMatthew Knepley   if (!nullspace) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
52649ef022SMatthew Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
53649ef022SMatthew Knepley   ierr = PetscDSSetObjective(ds, pfield, pressure_Private);CHKERRQ(ierr);
54649ef022SMatthew Knepley   ierr = MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs);CHKERRQ(ierr);
55649ef022SMatthew Knepley   if (Nv != 1) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %D", Nv);
56649ef022SMatthew Knepley   ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
57649ef022SMatthew Knepley   if (PetscAbsScalar(pintd) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
58649ef022SMatthew Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
59649ef022SMatthew Knepley   ierr = PetscMalloc2(Nf, &intc, Nf, &intn);CHKERRQ(ierr);
60649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx);CHKERRQ(ierr);
61649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
62649ef022SMatthew Knepley   ierr = VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]);CHKERRQ(ierr);
63649ef022SMatthew Knepley #if defined (PETSC_USE_DEBUG)
64649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
65649ef022SMatthew 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]));
66649ef022SMatthew Knepley #endif
67649ef022SMatthew Knepley   ierr = PetscFree2(intc, intn);CHKERRQ(ierr);
68649ef022SMatthew Knepley   PetscFunctionReturn(0);
69649ef022SMatthew Knepley }
70649ef022SMatthew Knepley 
71649ef022SMatthew Knepley /*@C
72649ef022SMatthew 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().
73649ef022SMatthew Knepley 
74649ef022SMatthew Knepley    Logically Collective on SNES
75649ef022SMatthew Knepley 
76649ef022SMatthew Knepley    Input Parameters:
77649ef022SMatthew Knepley +  snes - the SNES context
78649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
79649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
80649ef022SMatthew Knepley .  snorm - 2-norm of current step
81649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
82649ef022SMatthew Knepley -  ctx   - Optional user context
83649ef022SMatthew Knepley 
84649ef022SMatthew Knepley    Output Parameter:
85649ef022SMatthew Knepley .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
86649ef022SMatthew Knepley 
87649ef022SMatthew Knepley    Notes:
88649ef022SMatthew 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.
89649ef022SMatthew Knepley 
90649ef022SMatthew Knepley    Level: advanced
91649ef022SMatthew Knepley 
92649ef022SMatthew Knepley .seealso: SNESConvergedDefault(), SNESSetConvergenceTest(), DMSetNullSpaceConstructor()
93649ef022SMatthew Knepley @*/
94649ef022SMatthew Knepley PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
95649ef022SMatthew Knepley {
96649ef022SMatthew Knepley   PetscBool      monitorIntegral = PETSC_FALSE;
97649ef022SMatthew Knepley   PetscErrorCode ierr;
98649ef022SMatthew Knepley 
99649ef022SMatthew Knepley   PetscFunctionBegin;
100649ef022SMatthew Knepley   ierr = SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx);CHKERRQ(ierr);
101649ef022SMatthew Knepley   if (monitorIntegral) {
102649ef022SMatthew Knepley     Mat          J;
103649ef022SMatthew Knepley     Vec          u;
104649ef022SMatthew Knepley     MatNullSpace nullspace;
105649ef022SMatthew Knepley     const Vec   *nullvecs;
106649ef022SMatthew Knepley     PetscScalar  pintd;
107649ef022SMatthew Knepley 
108649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
109649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
110649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
111649ef022SMatthew Knepley     ierr = MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);CHKERRQ(ierr);
112649ef022SMatthew Knepley     ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
113649ef022SMatthew Knepley     ierr = PetscInfo1(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));CHKERRQ(ierr);
114649ef022SMatthew Knepley   }
115649ef022SMatthew Knepley   if (*reason > 0) {
116649ef022SMatthew Knepley     Mat          J;
117649ef022SMatthew Knepley     Vec          u;
118649ef022SMatthew Knepley     MatNullSpace nullspace;
119649ef022SMatthew Knepley     PetscInt     pfield = 1;
120649ef022SMatthew Knepley 
121649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
122649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
123649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
124649ef022SMatthew Knepley     ierr = SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx);CHKERRQ(ierr);
125649ef022SMatthew Knepley   }
126649ef022SMatthew Knepley   PetscFunctionReturn(0);
127649ef022SMatthew Knepley }
128649ef022SMatthew Knepley 
12924cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
13024cdb843SMatthew G. Knepley 
1316da023fcSToby Isaac static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
1326da023fcSToby Isaac {
1336da023fcSToby Isaac   PetscBool      isPlex;
1346da023fcSToby Isaac   PetscErrorCode ierr;
1356da023fcSToby Isaac 
1366da023fcSToby Isaac   PetscFunctionBegin;
1376da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
1386da023fcSToby Isaac   if (isPlex) {
1396da023fcSToby Isaac     *plex = dm;
1406da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
141f7148743SMatthew G. Knepley   } else {
142f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
143f7148743SMatthew G. Knepley     if (!*plex) {
1446da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
145f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
1466da023fcSToby Isaac       if (copy) {
1476da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
1489a2a23afSMatthew G. Knepley         ierr = DMCopyAuxiliaryVec(dm, *plex);CHKERRQ(ierr);
1496da023fcSToby Isaac       }
150f7148743SMatthew G. Knepley     } else {
151f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
152f7148743SMatthew G. Knepley     }
1536da023fcSToby Isaac   }
1546da023fcSToby Isaac   PetscFunctionReturn(0);
1556da023fcSToby Isaac }
1566da023fcSToby Isaac 
1574267b1a3SMatthew G. Knepley /*@C
1584267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
1594267b1a3SMatthew G. Knepley 
160d083f849SBarry Smith   Collective
1614267b1a3SMatthew G. Knepley 
1624267b1a3SMatthew G. Knepley   Input Parameter:
1634267b1a3SMatthew G. Knepley . comm - the communicator
1644267b1a3SMatthew G. Knepley 
1654267b1a3SMatthew G. Knepley   Output Parameter:
1664267b1a3SMatthew G. Knepley . ctx - the context
1674267b1a3SMatthew G. Knepley 
1684267b1a3SMatthew G. Knepley   Level: beginner
1694267b1a3SMatthew G. Knepley 
1704267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
1714267b1a3SMatthew G. Knepley @*/
1720adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
1730adebc6cSBarry Smith {
174552f7358SJed Brown   PetscErrorCode ierr;
175552f7358SJed Brown 
176552f7358SJed Brown   PetscFunctionBegin;
177552f7358SJed Brown   PetscValidPointer(ctx, 2);
17895dccacaSBarry Smith   ierr = PetscNew(ctx);CHKERRQ(ierr);
1791aa26658SKarl Rupp 
180552f7358SJed Brown   (*ctx)->comm   = comm;
181552f7358SJed Brown   (*ctx)->dim    = -1;
182552f7358SJed Brown   (*ctx)->nInput = 0;
1830298fd71SBarry Smith   (*ctx)->points = NULL;
1840298fd71SBarry Smith   (*ctx)->cells  = NULL;
185552f7358SJed Brown   (*ctx)->n      = -1;
1860298fd71SBarry Smith   (*ctx)->coords = NULL;
187552f7358SJed Brown   PetscFunctionReturn(0);
188552f7358SJed Brown }
189552f7358SJed Brown 
1904267b1a3SMatthew G. Knepley /*@C
1914267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
1924267b1a3SMatthew G. Knepley 
1934267b1a3SMatthew G. Knepley   Not collective
1944267b1a3SMatthew G. Knepley 
1954267b1a3SMatthew G. Knepley   Input Parameters:
1964267b1a3SMatthew G. Knepley + ctx - the context
1974267b1a3SMatthew G. Knepley - dim - the spatial dimension
1984267b1a3SMatthew G. Knepley 
1994267b1a3SMatthew G. Knepley   Level: intermediate
2004267b1a3SMatthew G. Knepley 
2014267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2024267b1a3SMatthew G. Knepley @*/
2030adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
2040adebc6cSBarry Smith {
205552f7358SJed Brown   PetscFunctionBegin;
206ff1e0c32SBarry Smith   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
207552f7358SJed Brown   ctx->dim = dim;
208552f7358SJed Brown   PetscFunctionReturn(0);
209552f7358SJed Brown }
210552f7358SJed Brown 
2114267b1a3SMatthew G. Knepley /*@C
2124267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2134267b1a3SMatthew G. Knepley 
2144267b1a3SMatthew G. Knepley   Not collective
2154267b1a3SMatthew G. Knepley 
2164267b1a3SMatthew G. Knepley   Input Parameter:
2174267b1a3SMatthew G. Knepley . ctx - the context
2184267b1a3SMatthew G. Knepley 
2194267b1a3SMatthew G. Knepley   Output Parameter:
2204267b1a3SMatthew G. Knepley . dim - the spatial dimension
2214267b1a3SMatthew G. Knepley 
2224267b1a3SMatthew G. Knepley   Level: intermediate
2234267b1a3SMatthew G. Knepley 
2244267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2254267b1a3SMatthew G. Knepley @*/
2260adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
2270adebc6cSBarry Smith {
228552f7358SJed Brown   PetscFunctionBegin;
229552f7358SJed Brown   PetscValidIntPointer(dim, 2);
230552f7358SJed Brown   *dim = ctx->dim;
231552f7358SJed Brown   PetscFunctionReturn(0);
232552f7358SJed Brown }
233552f7358SJed Brown 
2344267b1a3SMatthew G. Knepley /*@C
2354267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2364267b1a3SMatthew G. Knepley 
2374267b1a3SMatthew G. Knepley   Not collective
2384267b1a3SMatthew G. Knepley 
2394267b1a3SMatthew G. Knepley   Input Parameters:
2404267b1a3SMatthew G. Knepley + ctx - the context
2414267b1a3SMatthew G. Knepley - dof - the number of fields
2424267b1a3SMatthew G. Knepley 
2434267b1a3SMatthew G. Knepley   Level: intermediate
2444267b1a3SMatthew G. Knepley 
2454267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2464267b1a3SMatthew G. Knepley @*/
2470adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
2480adebc6cSBarry Smith {
249552f7358SJed Brown   PetscFunctionBegin;
250ff1e0c32SBarry Smith   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
251552f7358SJed Brown   ctx->dof = dof;
252552f7358SJed Brown   PetscFunctionReturn(0);
253552f7358SJed Brown }
254552f7358SJed Brown 
2554267b1a3SMatthew G. Knepley /*@C
2564267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2574267b1a3SMatthew G. Knepley 
2584267b1a3SMatthew G. Knepley   Not collective
2594267b1a3SMatthew G. Knepley 
2604267b1a3SMatthew G. Knepley   Input Parameter:
2614267b1a3SMatthew G. Knepley . ctx - the context
2624267b1a3SMatthew G. Knepley 
2634267b1a3SMatthew G. Knepley   Output Parameter:
2644267b1a3SMatthew G. Knepley . dof - the number of fields
2654267b1a3SMatthew G. Knepley 
2664267b1a3SMatthew G. Knepley   Level: intermediate
2674267b1a3SMatthew G. Knepley 
2684267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2694267b1a3SMatthew G. Knepley @*/
2700adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
2710adebc6cSBarry Smith {
272552f7358SJed Brown   PetscFunctionBegin;
273552f7358SJed Brown   PetscValidIntPointer(dof, 2);
274552f7358SJed Brown   *dof = ctx->dof;
275552f7358SJed Brown   PetscFunctionReturn(0);
276552f7358SJed Brown }
277552f7358SJed Brown 
2784267b1a3SMatthew G. Knepley /*@C
2794267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2804267b1a3SMatthew G. Knepley 
2814267b1a3SMatthew G. Knepley   Not collective
2824267b1a3SMatthew G. Knepley 
2834267b1a3SMatthew G. Knepley   Input Parameters:
2844267b1a3SMatthew G. Knepley + ctx    - the context
2854267b1a3SMatthew G. Knepley . n      - the number of points
2864267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2874267b1a3SMatthew G. Knepley 
2884267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
2894267b1a3SMatthew G. Knepley 
2904267b1a3SMatthew G. Knepley   Level: intermediate
2914267b1a3SMatthew G. Knepley 
2924267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
2934267b1a3SMatthew G. Knepley @*/
2940adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
2950adebc6cSBarry Smith {
296552f7358SJed Brown   PetscErrorCode ierr;
297552f7358SJed Brown 
298552f7358SJed Brown   PetscFunctionBegin;
2990adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
3000adebc6cSBarry Smith   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
301552f7358SJed Brown   ctx->nInput = n;
3021aa26658SKarl Rupp 
303785e854fSJed Brown   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
304580bdb30SBarry Smith   ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr);
305552f7358SJed Brown   PetscFunctionReturn(0);
306552f7358SJed Brown }
307552f7358SJed Brown 
3084267b1a3SMatthew G. Knepley /*@C
30952aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
3104267b1a3SMatthew G. Knepley 
3114267b1a3SMatthew G. Knepley   Collective on ctx
3124267b1a3SMatthew G. Knepley 
3134267b1a3SMatthew G. Knepley   Input Parameters:
3144267b1a3SMatthew G. Knepley + ctx - the context
3154267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
31652aa1562SMatthew G. Knepley . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
3177deeda50SSatish Balay - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
3184267b1a3SMatthew G. Knepley 
3194267b1a3SMatthew G. Knepley   Level: intermediate
3204267b1a3SMatthew G. Knepley 
3214267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3224267b1a3SMatthew G. Knepley @*/
32352aa1562SMatthew G. Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
3240adebc6cSBarry Smith {
325552f7358SJed Brown   MPI_Comm          comm = ctx->comm;
326552f7358SJed Brown   PetscScalar       *a;
327552f7358SJed Brown   PetscInt          p, q, i;
328552f7358SJed Brown   PetscMPIInt       rank, size;
329552f7358SJed Brown   PetscErrorCode    ierr;
330552f7358SJed Brown   Vec               pointVec;
3313a93e3b7SToby Isaac   PetscSF           cellSF;
332552f7358SJed Brown   PetscLayout       layout;
333552f7358SJed Brown   PetscReal         *globalPoints;
334cb313848SJed Brown   PetscScalar       *globalPointsScalar;
335552f7358SJed Brown   const PetscInt    *ranges;
336552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3373a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3383a93e3b7SToby Isaac   const PetscInt    *foundPoints;
339552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3403a93e3b7SToby Isaac   PetscInt          n, N, numFound;
341552f7358SJed Brown 
34219436ca2SJed Brown   PetscFunctionBegin;
343064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
344ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
345ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
3460adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
34719436ca2SJed Brown   /* Locate points */
34819436ca2SJed Brown   n = ctx->nInput;
349552f7358SJed Brown   if (!redundantPoints) {
350552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
351552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
352552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
353552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
354552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
355552f7358SJed Brown     /* Communicate all points to all processes */
356dcca6d9dSJed Brown     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
357552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
358552f7358SJed Brown     for (p = 0; p < size; ++p) {
359552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
360552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
361552f7358SJed Brown     }
362ffc4695bSBarry Smith     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRMPI(ierr);
363552f7358SJed Brown   } else {
364552f7358SJed Brown     N = n;
365552f7358SJed Brown     globalPoints = ctx->points;
36638ea73c8SJed Brown     counts = displs = NULL;
36738ea73c8SJed Brown     layout = NULL;
368552f7358SJed Brown   }
369552f7358SJed Brown #if 0
370dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
37119436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
372552f7358SJed Brown #else
373cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
37446b3086cSToby Isaac   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
37546b3086cSToby Isaac   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
376cb313848SJed Brown #else
377cb313848SJed Brown   globalPointsScalar = globalPoints;
378cb313848SJed Brown #endif
37904706141SMatthew G Knepley   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
380dcca6d9dSJed Brown   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
3817b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
3823a93e3b7SToby Isaac   cellSF = NULL;
3832d1fa6caSMatthew G. Knepley   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
3843a93e3b7SToby Isaac   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
385552f7358SJed Brown #endif
3863a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3873a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
388552f7358SJed Brown   }
389552f7358SJed Brown   /* Let the lowest rank process own each point */
390820f2d46SBarry Smith   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRMPI(ierr);
391552f7358SJed Brown   ctx->n = 0;
392552f7358SJed Brown   for (p = 0; p < N; ++p) {
39352aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
39452aa1562SMatthew 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));
39552aa1562SMatthew G. Knepley       else if (!rank) ++ctx->n;
39652aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
397552f7358SJed Brown   }
398552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
399785e854fSJed Brown   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
400552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
401552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
402552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
403c0dedaeaSBarry Smith   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
404552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
405552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
406552f7358SJed Brown     if (globalProcs[p] == rank) {
407552f7358SJed Brown       PetscInt d;
408552f7358SJed Brown 
4091aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
410f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
411f317b747SMatthew G. Knepley       ++q;
412552f7358SJed Brown     }
41352aa1562SMatthew G. Knepley     if (globalProcs[p] == size && !rank) {
41452aa1562SMatthew G. Knepley       PetscInt d;
41552aa1562SMatthew G. Knepley 
41652aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
41752aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
41852aa1562SMatthew G. Knepley       ++q;
41952aa1562SMatthew G. Knepley     }
420552f7358SJed Brown   }
421552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
422552f7358SJed Brown #if 0
423552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
424552f7358SJed Brown #else
425552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
4263a93e3b7SToby Isaac   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
427552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
428552f7358SJed Brown #endif
429cb313848SJed Brown   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
430d343d804SMatthew G. Knepley   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
431552f7358SJed Brown   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
432552f7358SJed Brown   PetscFunctionReturn(0);
433552f7358SJed Brown }
434552f7358SJed Brown 
4354267b1a3SMatthew G. Knepley /*@C
4364267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
4374267b1a3SMatthew G. Knepley 
4384267b1a3SMatthew G. Knepley   Collective on ctx
4394267b1a3SMatthew G. Knepley 
4404267b1a3SMatthew G. Knepley   Input Parameter:
4414267b1a3SMatthew G. Knepley . ctx - the context
4424267b1a3SMatthew G. Knepley 
4434267b1a3SMatthew G. Knepley   Output Parameter:
4444267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4454267b1a3SMatthew G. Knepley 
4464267b1a3SMatthew 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.
4474267b1a3SMatthew G. Knepley 
4484267b1a3SMatthew G. Knepley   Level: intermediate
4494267b1a3SMatthew G. Knepley 
4504267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4514267b1a3SMatthew G. Knepley @*/
4520adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
4530adebc6cSBarry Smith {
454552f7358SJed Brown   PetscFunctionBegin;
455552f7358SJed Brown   PetscValidPointer(coordinates, 2);
4560adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
457552f7358SJed Brown   *coordinates = ctx->coords;
458552f7358SJed Brown   PetscFunctionReturn(0);
459552f7358SJed Brown }
460552f7358SJed Brown 
4614267b1a3SMatthew G. Knepley /*@C
4624267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
4634267b1a3SMatthew G. Knepley 
4644267b1a3SMatthew G. Knepley   Collective on ctx
4654267b1a3SMatthew G. Knepley 
4664267b1a3SMatthew G. Knepley   Input Parameter:
4674267b1a3SMatthew G. Knepley . ctx - the context
4684267b1a3SMatthew G. Knepley 
4694267b1a3SMatthew G. Knepley   Output Parameter:
4704267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4714267b1a3SMatthew G. Knepley 
4724267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
4734267b1a3SMatthew G. Knepley 
4744267b1a3SMatthew G. Knepley   Level: intermediate
4754267b1a3SMatthew G. Knepley 
4764267b1a3SMatthew G. Knepley .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4774267b1a3SMatthew G. Knepley @*/
4780adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
4790adebc6cSBarry Smith {
480552f7358SJed Brown   PetscErrorCode ierr;
481552f7358SJed Brown 
482552f7358SJed Brown   PetscFunctionBegin;
483552f7358SJed Brown   PetscValidPointer(v, 2);
4840adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
485552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
486552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
487552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
488c0dedaeaSBarry Smith   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
489552f7358SJed Brown   PetscFunctionReturn(0);
490552f7358SJed Brown }
491552f7358SJed Brown 
4924267b1a3SMatthew G. Knepley /*@C
4934267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
4944267b1a3SMatthew G. Knepley 
4954267b1a3SMatthew G. Knepley   Collective on ctx
4964267b1a3SMatthew G. Knepley 
4974267b1a3SMatthew G. Knepley   Input Parameters:
4984267b1a3SMatthew G. Knepley + ctx - the context
4994267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
5004267b1a3SMatthew G. Knepley 
5014267b1a3SMatthew G. Knepley   Level: intermediate
5024267b1a3SMatthew G. Knepley 
5034267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
5044267b1a3SMatthew G. Knepley @*/
5050adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
5060adebc6cSBarry Smith {
507552f7358SJed Brown   PetscErrorCode ierr;
508552f7358SJed Brown 
509552f7358SJed Brown   PetscFunctionBegin;
510552f7358SJed Brown   PetscValidPointer(v, 2);
5110adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
512552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
513552f7358SJed Brown   PetscFunctionReturn(0);
514552f7358SJed Brown }
515552f7358SJed Brown 
5167a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
517a6dfd86eSKarl Rupp {
518552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
51956044e6dSMatthew G. Knepley   const PetscScalar *coords;
52056044e6dSMatthew G. Knepley   PetscScalar    *a;
521552f7358SJed Brown   PetscInt       p;
522552f7358SJed Brown   PetscErrorCode ierr;
523552f7358SJed Brown 
524552f7358SJed Brown   PetscFunctionBegin;
525dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
52656044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
527552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
528552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
529552f7358SJed Brown     PetscInt     c = ctx->cells[p];
530a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
531552f7358SJed Brown     PetscReal    xi[4];
532552f7358SJed Brown     PetscInt     d, f, comp;
533552f7358SJed Brown 
5348e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
535ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5360298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5371aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5381aa26658SKarl Rupp 
539552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
540552f7358SJed Brown       xi[d] = 0.0;
5411aa26658SKarl 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]);
5421aa26658SKarl 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];
543552f7358SJed Brown     }
5440298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
545552f7358SJed Brown   }
546552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
54756044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
548552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
549552f7358SJed Brown   PetscFunctionReturn(0);
550552f7358SJed Brown }
551552f7358SJed Brown 
5527a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
5537a1931ceSMatthew G. Knepley {
5547a1931ceSMatthew G. Knepley   PetscReal      *v0, *J, *invJ, detJ;
55556044e6dSMatthew G. Knepley   const PetscScalar *coords;
55656044e6dSMatthew G. Knepley   PetscScalar    *a;
5577a1931ceSMatthew G. Knepley   PetscInt       p;
5587a1931ceSMatthew G. Knepley   PetscErrorCode ierr;
5597a1931ceSMatthew G. Knepley 
5607a1931ceSMatthew G. Knepley   PetscFunctionBegin;
561dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
56256044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5637a1931ceSMatthew G. Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
5647a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5657a1931ceSMatthew G. Knepley     PetscInt       c = ctx->cells[p];
5667a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5672584bbe8SMatthew G. Knepley     PetscScalar   *x = NULL;
5687a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5697a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5707a1931ceSMatthew G. Knepley 
5718e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
572ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5737a1931ceSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5747a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5757a1931ceSMatthew G. Knepley 
5767a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5777a1931ceSMatthew G. Knepley       xi[d] = 0.0;
5787a1931ceSMatthew 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]);
5797a1931ceSMatthew 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];
5807a1931ceSMatthew G. Knepley     }
5817a1931ceSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5827a1931ceSMatthew G. Knepley   }
5837a1931ceSMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
58456044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5857a1931ceSMatthew G. Knepley   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
5867a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
5877a1931ceSMatthew G. Knepley }
5887a1931ceSMatthew G. Knepley 
5895820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
590552f7358SJed Brown {
591552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
592552f7358SJed Brown   const PetscScalar x0        = vertices[0];
593552f7358SJed Brown   const PetscScalar y0        = vertices[1];
594552f7358SJed Brown   const PetscScalar x1        = vertices[2];
595552f7358SJed Brown   const PetscScalar y1        = vertices[3];
596552f7358SJed Brown   const PetscScalar x2        = vertices[4];
597552f7358SJed Brown   const PetscScalar y2        = vertices[5];
598552f7358SJed Brown   const PetscScalar x3        = vertices[6];
599552f7358SJed Brown   const PetscScalar y3        = vertices[7];
600552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
601552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
602552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
603552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
604552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
605552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
60656044e6dSMatthew G. Knepley   const PetscScalar *ref;
60756044e6dSMatthew G. Knepley   PetscScalar       *real;
608552f7358SJed Brown   PetscErrorCode    ierr;
609552f7358SJed Brown 
610552f7358SJed Brown   PetscFunctionBegin;
61156044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
612552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
613552f7358SJed Brown   {
614552f7358SJed Brown     const PetscScalar p0 = ref[0];
615552f7358SJed Brown     const PetscScalar p1 = ref[1];
616552f7358SJed Brown 
617552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
618552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
619552f7358SJed Brown   }
620552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
62156044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
622552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
623552f7358SJed Brown   PetscFunctionReturn(0);
624552f7358SJed Brown }
625552f7358SJed Brown 
626af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
627d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
628552f7358SJed Brown {
629552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
630552f7358SJed Brown   const PetscScalar x0        = vertices[0];
631552f7358SJed Brown   const PetscScalar y0        = vertices[1];
632552f7358SJed Brown   const PetscScalar x1        = vertices[2];
633552f7358SJed Brown   const PetscScalar y1        = vertices[3];
634552f7358SJed Brown   const PetscScalar x2        = vertices[4];
635552f7358SJed Brown   const PetscScalar y2        = vertices[5];
636552f7358SJed Brown   const PetscScalar x3        = vertices[6];
637552f7358SJed Brown   const PetscScalar y3        = vertices[7];
638552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
639552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
64056044e6dSMatthew G. Knepley   const PetscScalar *ref;
641552f7358SJed Brown   PetscErrorCode    ierr;
642552f7358SJed Brown 
643552f7358SJed Brown   PetscFunctionBegin;
64456044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
645552f7358SJed Brown   {
646552f7358SJed Brown     const PetscScalar x       = ref[0];
647552f7358SJed Brown     const PetscScalar y       = ref[1];
648552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
649da80777bSKarl Rupp     PetscScalar       values[4];
650da80777bSKarl Rupp 
651da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
652da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
65394ab13aaSBarry Smith     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
654552f7358SJed Brown   }
655552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
65656044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
65794ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65894ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659552f7358SJed Brown   PetscFunctionReturn(0);
660552f7358SJed Brown }
661552f7358SJed Brown 
662a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
663a6dfd86eSKarl Rupp {
664fafc0619SMatthew G Knepley   DM             dmCoord;
665de73a395SMatthew G. Knepley   PetscFE        fem = NULL;
666552f7358SJed Brown   SNES           snes;
667552f7358SJed Brown   KSP            ksp;
668552f7358SJed Brown   PC             pc;
669552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
670552f7358SJed Brown   Mat            J;
671ef0bb6c7SMatthew G. Knepley   PetscTabulation    T;
67256044e6dSMatthew G. Knepley   const PetscScalar *coords;
67356044e6dSMatthew G. Knepley   PetscScalar    *a;
674ef0bb6c7SMatthew G. Knepley   PetscReal      xir[2];
675de73a395SMatthew G. Knepley   PetscInt       Nf, p;
6765509d985SMatthew G. Knepley   const PetscInt dof = ctx->dof;
677552f7358SJed Brown   PetscErrorCode ierr;
678552f7358SJed Brown 
679552f7358SJed Brown   PetscFunctionBegin;
680de73a395SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
68144a7f3ddSMatthew G. Knepley   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
682552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
683fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
684552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
685552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
686552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
687552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
688c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
689552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
690552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
691552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
692552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
693552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
694552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
6950298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
6960298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
697552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
698552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
699552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
700552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
701552f7358SJed Brown 
70256044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
703552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
704ef0bb6c7SMatthew G. Knepley   ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr);
705552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
706a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
707552f7358SJed Brown     PetscScalar *xi;
708552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
709552f7358SJed Brown 
710552f7358SJed Brown     /* Can make this do all points at once */
7110298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
712ff1e0c32SBarry Smith     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
7130298fd71SBarry Smith     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
7140298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
7150298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
716552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
717552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
718552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
719552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
720552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
721552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
722cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
723cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
7245509d985SMatthew G. Knepley     if (4*dof != xSize) {
7255509d985SMatthew G. Knepley       PetscInt d;
7261aa26658SKarl Rupp 
7275509d985SMatthew G. Knepley       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
728ef0bb6c7SMatthew G. Knepley       ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr);
7295509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7305509d985SMatthew G. Knepley         a[p*dof+comp] = 0.0;
7315509d985SMatthew G. Knepley         for (d = 0; d < xSize/dof; ++d) {
732ef0bb6c7SMatthew G. Knepley           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
7335509d985SMatthew G. Knepley         }
7345509d985SMatthew G. Knepley       }
7355509d985SMatthew G. Knepley     } else {
7365509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp)
7375509d985SMatthew 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];
7385509d985SMatthew G. Knepley     }
739552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
7400298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
7410298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
742552f7358SJed Brown   }
743ef0bb6c7SMatthew G. Knepley   ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
744552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
74556044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
746552f7358SJed Brown 
747552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
748552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
749552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
750552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
751552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
752552f7358SJed Brown   PetscFunctionReturn(0);
753552f7358SJed Brown }
754552f7358SJed Brown 
7555820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
756552f7358SJed Brown {
757552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
758552f7358SJed Brown   const PetscScalar x0        = vertices[0];
759552f7358SJed Brown   const PetscScalar y0        = vertices[1];
760552f7358SJed Brown   const PetscScalar z0        = vertices[2];
7617a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
7627a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
7637a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
764552f7358SJed Brown   const PetscScalar x2        = vertices[6];
765552f7358SJed Brown   const PetscScalar y2        = vertices[7];
766552f7358SJed Brown   const PetscScalar z2        = vertices[8];
7677a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
7687a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
7697a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
770552f7358SJed Brown   const PetscScalar x4        = vertices[12];
771552f7358SJed Brown   const PetscScalar y4        = vertices[13];
772552f7358SJed Brown   const PetscScalar z4        = vertices[14];
773552f7358SJed Brown   const PetscScalar x5        = vertices[15];
774552f7358SJed Brown   const PetscScalar y5        = vertices[16];
775552f7358SJed Brown   const PetscScalar z5        = vertices[17];
776552f7358SJed Brown   const PetscScalar x6        = vertices[18];
777552f7358SJed Brown   const PetscScalar y6        = vertices[19];
778552f7358SJed Brown   const PetscScalar z6        = vertices[20];
779552f7358SJed Brown   const PetscScalar x7        = vertices[21];
780552f7358SJed Brown   const PetscScalar y7        = vertices[22];
781552f7358SJed Brown   const PetscScalar z7        = vertices[23];
782552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
783552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
784552f7358SJed Brown   const PetscScalar h_1       = z1 - z0;
785552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
786552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
787552f7358SJed Brown   const PetscScalar h_3       = z3 - z0;
788552f7358SJed Brown   const PetscScalar f_4       = x4 - x0;
789552f7358SJed Brown   const PetscScalar g_4       = y4 - y0;
790552f7358SJed Brown   const PetscScalar h_4       = z4 - z0;
791552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
792552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
793552f7358SJed Brown   const PetscScalar h_01      = z2 - z1 - z3 + z0;
794552f7358SJed Brown   const PetscScalar f_12      = x7 - x3 - x4 + x0;
795552f7358SJed Brown   const PetscScalar g_12      = y7 - y3 - y4 + y0;
796552f7358SJed Brown   const PetscScalar h_12      = z7 - z3 - z4 + z0;
797552f7358SJed Brown   const PetscScalar f_02      = x5 - x1 - x4 + x0;
798552f7358SJed Brown   const PetscScalar g_02      = y5 - y1 - y4 + y0;
799552f7358SJed Brown   const PetscScalar h_02      = z5 - z1 - z4 + z0;
800552f7358SJed Brown   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
801552f7358SJed Brown   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
802552f7358SJed Brown   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
80356044e6dSMatthew G. Knepley   const PetscScalar *ref;
80456044e6dSMatthew G. Knepley   PetscScalar       *real;
805552f7358SJed Brown   PetscErrorCode    ierr;
806552f7358SJed Brown 
807552f7358SJed Brown   PetscFunctionBegin;
80856044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
809552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
810552f7358SJed Brown   {
811552f7358SJed Brown     const PetscScalar p0 = ref[0];
812552f7358SJed Brown     const PetscScalar p1 = ref[1];
813552f7358SJed Brown     const PetscScalar p2 = ref[2];
814552f7358SJed Brown 
815552f7358SJed 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;
816552f7358SJed 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;
817552f7358SJed 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;
818552f7358SJed Brown   }
819552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
82056044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
821552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
822552f7358SJed Brown   PetscFunctionReturn(0);
823552f7358SJed Brown }
824552f7358SJed Brown 
825d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
826552f7358SJed Brown {
827552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
828552f7358SJed Brown   const PetscScalar x0        = vertices[0];
829552f7358SJed Brown   const PetscScalar y0        = vertices[1];
830552f7358SJed Brown   const PetscScalar z0        = vertices[2];
8317a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
8327a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
8337a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
834552f7358SJed Brown   const PetscScalar x2        = vertices[6];
835552f7358SJed Brown   const PetscScalar y2        = vertices[7];
836552f7358SJed Brown   const PetscScalar z2        = vertices[8];
8377a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
8387a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
8397a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
840552f7358SJed Brown   const PetscScalar x4        = vertices[12];
841552f7358SJed Brown   const PetscScalar y4        = vertices[13];
842552f7358SJed Brown   const PetscScalar z4        = vertices[14];
843552f7358SJed Brown   const PetscScalar x5        = vertices[15];
844552f7358SJed Brown   const PetscScalar y5        = vertices[16];
845552f7358SJed Brown   const PetscScalar z5        = vertices[17];
846552f7358SJed Brown   const PetscScalar x6        = vertices[18];
847552f7358SJed Brown   const PetscScalar y6        = vertices[19];
848552f7358SJed Brown   const PetscScalar z6        = vertices[20];
849552f7358SJed Brown   const PetscScalar x7        = vertices[21];
850552f7358SJed Brown   const PetscScalar y7        = vertices[22];
851552f7358SJed Brown   const PetscScalar z7        = vertices[23];
852552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
853552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
854552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
855552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
856552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
857552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
858552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
859552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
860552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
861552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
862552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
863552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
86456044e6dSMatthew G. Knepley   const PetscScalar *ref;
865552f7358SJed Brown   PetscErrorCode    ierr;
866552f7358SJed Brown 
867552f7358SJed Brown   PetscFunctionBegin;
86856044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
869552f7358SJed Brown   {
870552f7358SJed Brown     const PetscScalar x       = ref[0];
871552f7358SJed Brown     const PetscScalar y       = ref[1];
872552f7358SJed Brown     const PetscScalar z       = ref[2];
873552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
874da80777bSKarl Rupp     PetscScalar       values[9];
875da80777bSKarl Rupp 
876da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
877da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
878da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
879da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
880da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
881da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
882da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
883da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
884da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
8851aa26658SKarl Rupp 
88694ab13aaSBarry Smith     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
887552f7358SJed Brown   }
888552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
88956044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
89094ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89194ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892552f7358SJed Brown   PetscFunctionReturn(0);
893552f7358SJed Brown }
894552f7358SJed Brown 
895a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
896a6dfd86eSKarl Rupp {
897fafc0619SMatthew G Knepley   DM             dmCoord;
898552f7358SJed Brown   SNES           snes;
899552f7358SJed Brown   KSP            ksp;
900552f7358SJed Brown   PC             pc;
901552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
902552f7358SJed Brown   Mat            J;
90356044e6dSMatthew G. Knepley   const PetscScalar *coords;
90456044e6dSMatthew G. Knepley   PetscScalar    *a;
905552f7358SJed Brown   PetscInt       p;
906552f7358SJed Brown   PetscErrorCode ierr;
907552f7358SJed Brown 
908552f7358SJed Brown   PetscFunctionBegin;
909552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
910fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
911552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
912552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
913552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
914552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
915c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
916552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
917552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
918552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
919552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
920552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
921552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
9220298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
9230298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
924552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
925552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
926552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
927552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
928552f7358SJed Brown 
92956044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
930552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
931552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
932a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
933552f7358SJed Brown     PetscScalar *xi;
934cb313848SJed Brown     PetscReal    xir[3];
935552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
936552f7358SJed Brown 
937552f7358SJed Brown     /* Can make this do all points at once */
9380298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
939ff1e0c32SBarry Smith     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
9400298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
941ff1e0c32SBarry Smith     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
9420298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
9430298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
944552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
945552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
946552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
947552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
948552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
949552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
950552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
951cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
952cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
953cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
954552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
955552f7358SJed Brown       a[p*ctx->dof+comp] =
956cb313848SJed Brown         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
9577a1931ceSMatthew G. Knepley         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
958cb313848SJed Brown         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
9597a1931ceSMatthew G. Knepley         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
960cb313848SJed Brown         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
961cb313848SJed Brown         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
962cb313848SJed Brown         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
963cb313848SJed Brown         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
964552f7358SJed Brown     }
965552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
9660298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
9670298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
968552f7358SJed Brown   }
969552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
97056044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
971552f7358SJed Brown 
972552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
973552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
974552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
975552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
976552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
977552f7358SJed Brown   PetscFunctionReturn(0);
978552f7358SJed Brown }
979552f7358SJed Brown 
9804267b1a3SMatthew G. Knepley /*@C
9814267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9824267b1a3SMatthew G. Knepley 
983552f7358SJed Brown   Input Parameters:
984552f7358SJed Brown + ctx - The DMInterpolationInfo context
985552f7358SJed Brown . dm  - The DM
986552f7358SJed Brown - x   - The local vector containing the field to be interpolated
987552f7358SJed Brown 
988552f7358SJed Brown   Output Parameters:
989552f7358SJed Brown . v   - The vector containing the interpolated values
9904267b1a3SMatthew G. Knepley 
9914267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
9924267b1a3SMatthew G. Knepley 
9934267b1a3SMatthew G. Knepley   Level: beginner
9944267b1a3SMatthew G. Knepley 
9954267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
9964267b1a3SMatthew G. Knepley @*/
9970adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
9980adebc6cSBarry Smith {
999680d18d5SMatthew G. Knepley   PetscInt       n;
1000552f7358SJed Brown   PetscErrorCode ierr;
1001552f7358SJed Brown 
1002552f7358SJed Brown   PetscFunctionBegin;
1003552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1004552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1005552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
1006552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1007ff1e0c32SBarry 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);
1008552f7358SJed Brown   if (n) {
1009680d18d5SMatthew G. Knepley     PetscDS        ds;
1010680d18d5SMatthew G. Knepley     DMPolytopeType ct;
1011680d18d5SMatthew G. Knepley     PetscBool      done = PETSC_FALSE;
1012680d18d5SMatthew G. Knepley 
1013680d18d5SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1014680d18d5SMatthew G. Knepley     if (ds) {
1015680d18d5SMatthew G. Knepley       const PetscScalar *coords;
1016680d18d5SMatthew G. Knepley       PetscScalar       *interpolant;
1017680d18d5SMatthew G. Knepley       PetscInt           cdim, d, p, Nf, field, c = 0;
1018680d18d5SMatthew G. Knepley 
1019680d18d5SMatthew G. Knepley       ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1020680d18d5SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
1021680d18d5SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
1022680d18d5SMatthew G. Knepley         PetscTabulation T;
1023680d18d5SMatthew G. Knepley         PetscFE         fe;
1024680d18d5SMatthew G. Knepley         PetscClassId    id;
1025680d18d5SMatthew G. Knepley         PetscReal       xi[3];
1026680d18d5SMatthew G. Knepley         PetscInt        Nc, f, fc;
1027680d18d5SMatthew G. Knepley 
1028680d18d5SMatthew G. Knepley         ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
1029680d18d5SMatthew G. Knepley         ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
1030680d18d5SMatthew G. Knepley         if (id != PETSCFE_CLASSID) break;
1031680d18d5SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1032680d18d5SMatthew G. Knepley         ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1033680d18d5SMatthew G. Knepley         ierr = VecGetArrayWrite(v, &interpolant);CHKERRQ(ierr);
1034680d18d5SMatthew G. Knepley         for (p = 0; p < ctx->n; ++p) {
1035680d18d5SMatthew G. Knepley           PetscScalar *xa = NULL;
1036680d18d5SMatthew G. Knepley           PetscReal    pcoords[3];
1037680d18d5SMatthew G. Knepley 
103852aa1562SMatthew G. Knepley           if (ctx->cells[p] < 0) continue;
1039680d18d5SMatthew G. Knepley           for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
1040680d18d5SMatthew G. Knepley           ierr = DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi);CHKERRQ(ierr);
1041680d18d5SMatthew G. Knepley           ierr = DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1042680d18d5SMatthew G. Knepley           ierr = PetscFECreateTabulation(fe, 1, 1, xi, 0, &T);CHKERRQ(ierr);
1043680d18d5SMatthew G. Knepley           {
1044680d18d5SMatthew G. Knepley             const PetscReal *basis = T->T[0];
1045680d18d5SMatthew G. Knepley             const PetscInt   Nb    = T->Nb;
1046680d18d5SMatthew G. Knepley             const PetscInt   Nc    = T->Nc;
1047680d18d5SMatthew G. Knepley             for (fc = 0; fc < Nc; ++fc) {
1048680d18d5SMatthew G. Knepley               interpolant[p*ctx->dof+c+fc] = 0.0;
1049680d18d5SMatthew G. Knepley               for (f = 0; f < Nb; ++f) {
1050680d18d5SMatthew G. Knepley                 interpolant[p*ctx->dof+c+fc] += xa[f]*basis[(0*Nb + f)*Nc + fc];
1051552f7358SJed Brown               }
1052680d18d5SMatthew G. Knepley             }
1053680d18d5SMatthew G. Knepley           }
1054680d18d5SMatthew G. Knepley           ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
1055680d18d5SMatthew G. Knepley           ierr = DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1056680d18d5SMatthew G. Knepley         }
1057680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayWrite(v, &interpolant);CHKERRQ(ierr);
1058680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1059680d18d5SMatthew G. Knepley         c += Nc;
1060680d18d5SMatthew G. Knepley       }
1061680d18d5SMatthew G. Knepley       if (field == Nf) {
1062680d18d5SMatthew G. Knepley         done = PETSC_TRUE;
1063680d18d5SMatthew G. Knepley         if (c != ctx->dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %D != %D dof specified for interpolation", c, ctx->dof);
1064680d18d5SMatthew G. Knepley       }
1065680d18d5SMatthew G. Knepley     }
1066680d18d5SMatthew G. Knepley     if (!done) {
1067680d18d5SMatthew G. Knepley       /* TODO Check each cell individually */
1068680d18d5SMatthew G. Knepley       ierr = DMPlexGetCellType(dm, ctx->cells[0], &ct);CHKERRQ(ierr);
1069680d18d5SMatthew G. Knepley       switch (ct) {
1070680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:      ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1071680d18d5SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL: ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1072680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TETRAHEDRON:   ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1073680d18d5SMatthew G. Knepley         case DM_POLYTOPE_HEXAHEDRON:    ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1074680d18d5SMatthew G. Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support fpr cell type %s", DMPolytopeTypes[ct]);
1075680d18d5SMatthew G. Knepley       }
1076680d18d5SMatthew G. Knepley     }
1077552f7358SJed Brown   }
1078552f7358SJed Brown   PetscFunctionReturn(0);
1079552f7358SJed Brown }
1080552f7358SJed Brown 
10814267b1a3SMatthew G. Knepley /*@C
10824267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
10834267b1a3SMatthew G. Knepley 
10844267b1a3SMatthew G. Knepley   Collective on ctx
10854267b1a3SMatthew G. Knepley 
10864267b1a3SMatthew G. Knepley   Input Parameter:
10874267b1a3SMatthew G. Knepley . ctx - the context
10884267b1a3SMatthew G. Knepley 
10894267b1a3SMatthew G. Knepley   Level: beginner
10904267b1a3SMatthew G. Knepley 
10914267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
10924267b1a3SMatthew G. Knepley @*/
10930adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
10940adebc6cSBarry Smith {
1095552f7358SJed Brown   PetscErrorCode ierr;
1096552f7358SJed Brown 
1097552f7358SJed Brown   PetscFunctionBegin;
1098064a246eSJacob Faibussowitsch   PetscValidPointer(ctx, 1);
1099552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
1100552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
1101552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
1102552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
11030298fd71SBarry Smith   *ctx = NULL;
1104552f7358SJed Brown   PetscFunctionReturn(0);
1105552f7358SJed Brown }
1106cc0c4584SMatthew G. Knepley 
1107cc0c4584SMatthew G. Knepley /*@C
1108cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1109cc0c4584SMatthew G. Knepley 
1110cc0c4584SMatthew G. Knepley   Collective on SNES
1111cc0c4584SMatthew G. Knepley 
1112cc0c4584SMatthew G. Knepley   Input Parameters:
1113cc0c4584SMatthew G. Knepley + snes   - the SNES context
1114cc0c4584SMatthew G. Knepley . its    - iteration number
1115cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1116d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
1117cc0c4584SMatthew G. Knepley 
1118cc0c4584SMatthew G. Knepley   Notes:
1119cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1120cc0c4584SMatthew G. Knepley 
1121cc0c4584SMatthew G. Knepley   Level: intermediate
1122cc0c4584SMatthew G. Knepley 
1123cc0c4584SMatthew G. Knepley .seealso: SNESMonitorSet(), SNESMonitorDefault()
1124cc0c4584SMatthew G. Knepley @*/
1125d43b4f6eSBarry Smith PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1126cc0c4584SMatthew G. Knepley {
1127d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1128cc0c4584SMatthew G. Knepley   Vec                res;
1129cc0c4584SMatthew G. Knepley   DM                 dm;
1130cc0c4584SMatthew G. Knepley   PetscSection       s;
1131cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1132cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1133cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1134cc0c4584SMatthew G. Knepley   PetscErrorCode     ierr;
1135cc0c4584SMatthew G. Knepley 
1136cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11374d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
11389e5d0892SLisandro Dalcin   ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr);
1139cc0c4584SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
114092fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
1141cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
1142cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1143cc0c4584SMatthew G. Knepley   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
1144cc0c4584SMatthew G. Knepley   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
1145cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1146cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1147cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1148cc0c4584SMatthew G. Knepley 
1149cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1150cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1151cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1152cc0c4584SMatthew G. Knepley     }
1153cc0c4584SMatthew G. Knepley   }
1154cc0c4584SMatthew G. Knepley   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
1155820f2d46SBarry Smith   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRMPI(ierr);
1156d43b4f6eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1157cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1158cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
1159cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
1160cc0c4584SMatthew G. Knepley     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
1161cc0c4584SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
1162cc0c4584SMatthew G. Knepley   }
1163cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
1164cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1165d43b4f6eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1166cc0c4584SMatthew G. Knepley   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
1167cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1168cc0c4584SMatthew G. Knepley }
116924cdb843SMatthew G. Knepley 
117024cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
117124cdb843SMatthew G. Knepley 
11726528b96dSMatthew G. Knepley PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
11736528b96dSMatthew G. Knepley {
11746528b96dSMatthew G. Knepley   PetscInt       depth;
11756528b96dSMatthew G. Knepley   PetscErrorCode ierr;
11766528b96dSMatthew G. Knepley 
11776528b96dSMatthew G. Knepley   PetscFunctionBegin;
11786528b96dSMatthew G. Knepley   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
11796528b96dSMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, cellIS);CHKERRQ(ierr);
11806528b96dSMatthew G. Knepley   if (!*cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, cellIS);CHKERRQ(ierr);}
11816528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
11826528b96dSMatthew G. Knepley }
11836528b96dSMatthew G. Knepley 
118424cdb843SMatthew G. Knepley /*@
118524cdb843SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
118624cdb843SMatthew G. Knepley 
118724cdb843SMatthew G. Knepley   Input Parameters:
118824cdb843SMatthew G. Knepley + dm - The mesh
118924cdb843SMatthew G. Knepley . X  - Local solution
119024cdb843SMatthew G. Knepley - user - The user context
119124cdb843SMatthew G. Knepley 
119224cdb843SMatthew G. Knepley   Output Parameter:
119324cdb843SMatthew G. Knepley . F  - Local output vector
119424cdb843SMatthew G. Knepley 
119524cdb843SMatthew G. Knepley   Level: developer
119624cdb843SMatthew G. Knepley 
11977a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
119824cdb843SMatthew G. Knepley @*/
119924cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
120024cdb843SMatthew G. Knepley {
12016da023fcSToby Isaac   DM             plex;
1202083401c6SMatthew G. Knepley   IS             allcellIS;
12036528b96dSMatthew G. Knepley   PetscInt       Nds, s;
120424cdb843SMatthew G. Knepley   PetscErrorCode ierr;
120524cdb843SMatthew G. Knepley 
120624cdb843SMatthew G. Knepley   PetscFunctionBegin;
12076da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
12086528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
12096528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
12106528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
12116528b96dSMatthew G. Knepley     PetscDS          ds;
12126528b96dSMatthew G. Knepley     IS               cellIS;
12136528b96dSMatthew G. Knepley     PetscHashFormKey key;
12146528b96dSMatthew G. Knepley 
12156528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
12166528b96dSMatthew G. Knepley     key.value = 0;
12176528b96dSMatthew G. Knepley     key.field = 0;
12186528b96dSMatthew G. Knepley     if (!key.label) {
12196528b96dSMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
12206528b96dSMatthew G. Knepley       cellIS = allcellIS;
12216528b96dSMatthew G. Knepley     } else {
12226528b96dSMatthew G. Knepley       IS pointIS;
12236528b96dSMatthew G. Knepley 
12246528b96dSMatthew G. Knepley       key.value = 1;
12256528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
12266528b96dSMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
12276528b96dSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
12286528b96dSMatthew G. Knepley     }
12296528b96dSMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
12306528b96dSMatthew G. Knepley     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
12316528b96dSMatthew G. Knepley   }
12326528b96dSMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
12336528b96dSMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
12346528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
12356528b96dSMatthew G. Knepley }
12366528b96dSMatthew G. Knepley 
12376528b96dSMatthew G. Knepley PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user)
12386528b96dSMatthew G. Knepley {
12396528b96dSMatthew G. Knepley   DM             plex;
12406528b96dSMatthew G. Knepley   IS             allcellIS;
12416528b96dSMatthew G. Knepley   PetscInt       Nds, s;
12426528b96dSMatthew G. Knepley   PetscErrorCode ierr;
12436528b96dSMatthew G. Knepley 
12446528b96dSMatthew G. Knepley   PetscFunctionBegin;
12456528b96dSMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
12466528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
12476528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1248083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1249083401c6SMatthew G. Knepley     PetscDS ds;
1250083401c6SMatthew G. Knepley     DMLabel label;
1251083401c6SMatthew G. Knepley     IS      cellIS;
1252083401c6SMatthew G. Knepley 
1253083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
12546528b96dSMatthew G. Knepley     {
12556528b96dSMatthew G. Knepley       PetscHMapForm     resmap[2] = {ds->wf->f0, ds->wf->f1};
12566528b96dSMatthew G. Knepley       PetscWeakForm     wf;
12576528b96dSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
12586528b96dSMatthew G. Knepley       PetscHashFormKey *reskeys;
12596528b96dSMatthew G. Knepley 
12606528b96dSMatthew G. Knepley       /* Get unique residual keys */
12616528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
12626528b96dSMatthew G. Knepley         PetscInt Nkm;
12636528b96dSMatthew G. Knepley         ierr = PetscHMapFormGetSize(resmap[m], &Nkm);CHKERRQ(ierr);
12646528b96dSMatthew G. Knepley         Nk  += Nkm;
12656528b96dSMatthew G. Knepley       }
12666528b96dSMatthew G. Knepley       ierr = PetscMalloc1(Nk, &reskeys);CHKERRQ(ierr);
12676528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
12686528b96dSMatthew G. Knepley         ierr = PetscHMapFormGetKeys(resmap[m], &off, reskeys);CHKERRQ(ierr);
12696528b96dSMatthew G. Knepley       }
12706528b96dSMatthew G. Knepley       if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
12716528b96dSMatthew G. Knepley       ierr = PetscHashFormKeySort(Nk, reskeys);CHKERRQ(ierr);
12726528b96dSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
12736528b96dSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
12746528b96dSMatthew G. Knepley           ++k;
12756528b96dSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
12766528b96dSMatthew G. Knepley         }
12776528b96dSMatthew G. Knepley       }
12786528b96dSMatthew G. Knepley       Nk = k;
12796528b96dSMatthew G. Knepley 
12806528b96dSMatthew G. Knepley       ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
12816528b96dSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
12826528b96dSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
12836528b96dSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
12846528b96dSMatthew G. Knepley 
1285083401c6SMatthew G. Knepley         if (!label) {
1286083401c6SMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1287083401c6SMatthew G. Knepley           cellIS = allcellIS;
1288083401c6SMatthew G. Knepley         } else {
1289083401c6SMatthew G. Knepley           IS pointIS;
1290083401c6SMatthew G. Knepley 
12916528b96dSMatthew G. Knepley           ierr = DMLabelGetStratumIS(label, val, &pointIS);CHKERRQ(ierr);
1292083401c6SMatthew G. Knepley           ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1293083401c6SMatthew G. Knepley           ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
12944a3e9fdbSToby Isaac         }
12956528b96dSMatthew G. Knepley         ierr = DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
12964a3e9fdbSToby Isaac         ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1297083401c6SMatthew G. Knepley       }
12986528b96dSMatthew G. Knepley       ierr = PetscFree(reskeys);CHKERRQ(ierr);
12996528b96dSMatthew G. Knepley     }
13006528b96dSMatthew G. Knepley   }
1301083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
13029a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
130324cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
130424cdb843SMatthew G. Knepley }
130524cdb843SMatthew G. Knepley 
1306bdd6f66aSToby Isaac /*@
1307bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1308bdd6f66aSToby Isaac 
1309bdd6f66aSToby Isaac   Input Parameters:
1310bdd6f66aSToby Isaac + dm - The mesh
1311bdd6f66aSToby Isaac - user - The user context
1312bdd6f66aSToby Isaac 
1313bdd6f66aSToby Isaac   Output Parameter:
1314bdd6f66aSToby Isaac . X  - Local solution
1315bdd6f66aSToby Isaac 
1316bdd6f66aSToby Isaac   Level: developer
1317bdd6f66aSToby Isaac 
13187a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
1319bdd6f66aSToby Isaac @*/
1320bdd6f66aSToby Isaac PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1321bdd6f66aSToby Isaac {
1322bdd6f66aSToby Isaac   DM             plex;
1323bdd6f66aSToby Isaac   PetscErrorCode ierr;
1324bdd6f66aSToby Isaac 
1325bdd6f66aSToby Isaac   PetscFunctionBegin;
1326bdd6f66aSToby Isaac   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1327bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1328bdd6f66aSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1329bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1330bdd6f66aSToby Isaac }
1331bdd6f66aSToby Isaac 
13327a73cf09SMatthew G. Knepley /*@
1333*8e3a2eefSMatthew G. Knepley   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
13347a73cf09SMatthew G. Knepley 
13357a73cf09SMatthew G. Knepley   Input Parameters:
1336*8e3a2eefSMatthew G. Knepley + dm   - The DM
13377a73cf09SMatthew G. Knepley . X    - Local solution vector
13387a73cf09SMatthew G. Knepley . Y    - Local input vector
13397a73cf09SMatthew G. Knepley - user - The user context
13407a73cf09SMatthew G. Knepley 
13417a73cf09SMatthew G. Knepley   Output Parameter:
1342*8e3a2eefSMatthew G. Knepley . F    - lcoal output vector
13437a73cf09SMatthew G. Knepley 
13447a73cf09SMatthew G. Knepley   Level: developer
13457a73cf09SMatthew G. Knepley 
1346*8e3a2eefSMatthew G. Knepley   Notes:
1347*8e3a2eefSMatthew G. Knepley   Users will typically use DMSNESCreateJacobianMF() followed by MatMult() instead of calling this routine directly.
1348*8e3a2eefSMatthew G. Knepley 
1349*8e3a2eefSMatthew G. Knepley .seealso: DMSNESCreateJacobianMF()
13507a73cf09SMatthew G. Knepley @*/
1351*8e3a2eefSMatthew G. Knepley PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
1352a925c78cSMatthew G. Knepley {
1353*8e3a2eefSMatthew G. Knepley   DM             plex;
1354*8e3a2eefSMatthew G. Knepley   IS             allcellIS;
1355*8e3a2eefSMatthew G. Knepley   PetscInt       Nds, s;
1356a925c78cSMatthew G. Knepley   PetscErrorCode ierr;
1357a925c78cSMatthew G. Knepley 
1358a925c78cSMatthew G. Knepley   PetscFunctionBegin;
13597a73cf09SMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1360*8e3a2eefSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1361*8e3a2eefSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1362*8e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1363*8e3a2eefSMatthew G. Knepley     PetscDS ds;
1364*8e3a2eefSMatthew G. Knepley     DMLabel label;
1365*8e3a2eefSMatthew G. Knepley     IS      cellIS;
13667a73cf09SMatthew G. Knepley 
1367*8e3a2eefSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1368*8e3a2eefSMatthew G. Knepley     {
1369*8e3a2eefSMatthew G. Knepley       PetscHMapForm     jacmap[4] = {ds->wf->g0, ds->wf->g1, ds->wf->g2, ds->wf->g3};
1370*8e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
1371*8e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
1372*8e3a2eefSMatthew G. Knepley       PetscHashFormKey *jackeys;
1373*8e3a2eefSMatthew G. Knepley 
1374*8e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
1375*8e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
1376*8e3a2eefSMatthew G. Knepley         PetscInt Nkm;
1377*8e3a2eefSMatthew G. Knepley         ierr = PetscHMapFormGetSize(jacmap[m], &Nkm);CHKERRQ(ierr);
1378*8e3a2eefSMatthew G. Knepley         Nk  += Nkm;
1379*8e3a2eefSMatthew G. Knepley       }
1380*8e3a2eefSMatthew G. Knepley       ierr = PetscMalloc1(Nk, &jackeys);CHKERRQ(ierr);
1381*8e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
1382*8e3a2eefSMatthew G. Knepley         ierr = PetscHMapFormGetKeys(jacmap[m], &off, jackeys);CHKERRQ(ierr);
1383*8e3a2eefSMatthew G. Knepley       }
1384*8e3a2eefSMatthew G. Knepley       if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1385*8e3a2eefSMatthew G. Knepley       ierr = PetscHashFormKeySort(Nk, jackeys);CHKERRQ(ierr);
1386*8e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
1387*8e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
1388*8e3a2eefSMatthew G. Knepley           ++k;
1389*8e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
1390*8e3a2eefSMatthew G. Knepley         }
1391*8e3a2eefSMatthew G. Knepley       }
1392*8e3a2eefSMatthew G. Knepley       Nk = k;
1393*8e3a2eefSMatthew G. Knepley 
1394*8e3a2eefSMatthew G. Knepley       ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1395*8e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
1396*8e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
1397*8e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
1398*8e3a2eefSMatthew G. Knepley 
1399*8e3a2eefSMatthew G. Knepley         if (!label) {
1400*8e3a2eefSMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1401*8e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
14027a73cf09SMatthew G. Knepley         } else {
1403*8e3a2eefSMatthew G. Knepley           IS pointIS;
1404a925c78cSMatthew G. Knepley 
1405*8e3a2eefSMatthew G. Knepley           ierr = DMLabelGetStratumIS(label, val, &pointIS);CHKERRQ(ierr);
1406*8e3a2eefSMatthew G. Knepley           ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1407*8e3a2eefSMatthew G. Knepley           ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1408a925c78cSMatthew G. Knepley         }
1409*8e3a2eefSMatthew G. Knepley         ierr = DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user);CHKERRQ(ierr);
14107a73cf09SMatthew G. Knepley         ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1411*8e3a2eefSMatthew G. Knepley       }
1412*8e3a2eefSMatthew G. Knepley       ierr = PetscFree(jackeys);CHKERRQ(ierr);
1413*8e3a2eefSMatthew G. Knepley     }
1414*8e3a2eefSMatthew G. Knepley   }
1415*8e3a2eefSMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
14167a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1417a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1418a925c78cSMatthew G. Knepley }
1419a925c78cSMatthew G. Knepley 
142024cdb843SMatthew G. Knepley /*@
142124cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
142224cdb843SMatthew G. Knepley 
142324cdb843SMatthew G. Knepley   Input Parameters:
142424cdb843SMatthew G. Knepley + dm - The mesh
142524cdb843SMatthew G. Knepley . X  - Local input vector
142624cdb843SMatthew G. Knepley - user - The user context
142724cdb843SMatthew G. Knepley 
142824cdb843SMatthew G. Knepley   Output Parameter:
142924cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
143024cdb843SMatthew G. Knepley 
143124cdb843SMatthew G. Knepley   Note:
143224cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
143324cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
143424cdb843SMatthew G. Knepley 
143524cdb843SMatthew G. Knepley   Level: developer
143624cdb843SMatthew G. Knepley 
143724cdb843SMatthew G. Knepley .seealso: FormFunctionLocal()
143824cdb843SMatthew G. Knepley @*/
143924cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
144024cdb843SMatthew G. Knepley {
14416da023fcSToby Isaac   DM             plex;
1442083401c6SMatthew G. Knepley   IS             allcellIS;
1443f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
14446528b96dSMatthew G. Knepley   PetscInt       Nds, s;
144524cdb843SMatthew G. Knepley   PetscErrorCode ierr;
144624cdb843SMatthew G. Knepley 
144724cdb843SMatthew G. Knepley   PetscFunctionBegin;
14486da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
14496528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
14506528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1451083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1452083401c6SMatthew G. Knepley     PetscDS          ds;
1453083401c6SMatthew G. Knepley     IS               cellIS;
14546528b96dSMatthew G. Knepley     PetscHashFormKey key;
1455083401c6SMatthew G. Knepley 
14566528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
14576528b96dSMatthew G. Knepley     key.value = 0;
14586528b96dSMatthew G. Knepley     key.field = 0;
14596528b96dSMatthew G. Knepley     if (!key.label) {
1460083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1461083401c6SMatthew G. Knepley       cellIS = allcellIS;
1462083401c6SMatthew G. Knepley     } else {
1463083401c6SMatthew G. Knepley       IS pointIS;
1464083401c6SMatthew G. Knepley 
14656528b96dSMatthew G. Knepley       key.value = 1;
14666528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1467083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1468083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1469083401c6SMatthew G. Knepley     }
1470083401c6SMatthew G. Knepley     if (!s) {
1471083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1472083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1473f04eb4edSMatthew G. Knepley       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
1474f04eb4edSMatthew G. Knepley       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1475083401c6SMatthew G. Knepley     }
14766528b96dSMatthew G. Knepley     ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
14774a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1478083401c6SMatthew G. Knepley   }
1479083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
14809a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
148124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
148224cdb843SMatthew G. Knepley }
14831878804aSMatthew G. Knepley 
1484*8e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx
1485*8e3a2eefSMatthew G. Knepley {
1486*8e3a2eefSMatthew G. Knepley   DM    dm;
1487*8e3a2eefSMatthew G. Knepley   Vec   X;
1488*8e3a2eefSMatthew G. Knepley   void *ctx;
1489*8e3a2eefSMatthew G. Knepley };
1490*8e3a2eefSMatthew G. Knepley 
1491*8e3a2eefSMatthew G. Knepley static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
1492*8e3a2eefSMatthew G. Knepley {
1493*8e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
1494*8e3a2eefSMatthew G. Knepley   PetscErrorCode               ierr;
1495*8e3a2eefSMatthew G. Knepley 
1496*8e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
1497*8e3a2eefSMatthew G. Knepley   ierr = MatShellGetContext(A, (void **) &ctx);CHKERRQ(ierr);
1498*8e3a2eefSMatthew G. Knepley   ierr = MatShellSetContext(A, NULL);CHKERRQ(ierr);
1499*8e3a2eefSMatthew G. Knepley   ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr);
1500*8e3a2eefSMatthew G. Knepley   ierr = VecDestroy(&ctx->X);CHKERRQ(ierr);
1501*8e3a2eefSMatthew G. Knepley   ierr = PetscFree(ctx);CHKERRQ(ierr);
1502*8e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
1503*8e3a2eefSMatthew G. Knepley }
1504*8e3a2eefSMatthew G. Knepley 
1505*8e3a2eefSMatthew G. Knepley static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
1506*8e3a2eefSMatthew G. Knepley {
1507*8e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
1508*8e3a2eefSMatthew G. Knepley   PetscErrorCode               ierr;
1509*8e3a2eefSMatthew G. Knepley 
1510*8e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
1511*8e3a2eefSMatthew G. Knepley   ierr = MatShellGetContext(A, (void **) &ctx);CHKERRQ(ierr);
1512*8e3a2eefSMatthew G. Knepley   ierr = DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx);CHKERRQ(ierr);
1513*8e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
1514*8e3a2eefSMatthew G. Knepley }
1515*8e3a2eefSMatthew G. Knepley 
1516*8e3a2eefSMatthew G. Knepley /*@
1517*8e3a2eefSMatthew G. Knepley   DMSNESCreateJacobianMF - Create a Mat which computes the action of the Jacobian matrix-free
1518*8e3a2eefSMatthew G. Knepley 
1519*8e3a2eefSMatthew G. Knepley   Collective on dm
1520*8e3a2eefSMatthew G. Knepley 
1521*8e3a2eefSMatthew G. Knepley   Input Parameters:
1522*8e3a2eefSMatthew G. Knepley + dm   - The DM
1523*8e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
1524*8e3a2eefSMatthew G. Knepley - user - A user context, or NULL
1525*8e3a2eefSMatthew G. Knepley 
1526*8e3a2eefSMatthew G. Knepley   Output Parameter:
1527*8e3a2eefSMatthew G. Knepley . J    - The Mat
1528*8e3a2eefSMatthew G. Knepley 
1529*8e3a2eefSMatthew G. Knepley   Level: advanced
1530*8e3a2eefSMatthew G. Knepley 
1531*8e3a2eefSMatthew G. Knepley   Notes:
1532*8e3a2eefSMatthew G. Knepley   Vec X is kept in Mat J, so updating X then updates the evaluation point.
1533*8e3a2eefSMatthew G. Knepley 
1534*8e3a2eefSMatthew G. Knepley .seealso: DMSNESComputeJacobianAction()
1535*8e3a2eefSMatthew G. Knepley @*/
1536*8e3a2eefSMatthew G. Knepley PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
1537*8e3a2eefSMatthew G. Knepley {
1538*8e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
1539*8e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
1540*8e3a2eefSMatthew G. Knepley   PetscErrorCode               ierr;
1541*8e3a2eefSMatthew G. Knepley 
1542*8e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
1543*8e3a2eefSMatthew G. Knepley   ierr = MatCreate(PetscObjectComm((PetscObject) dm), J);CHKERRQ(ierr);
1544*8e3a2eefSMatthew G. Knepley   ierr = MatSetType(*J, MATSHELL);CHKERRQ(ierr);
1545*8e3a2eefSMatthew G. Knepley   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
1546*8e3a2eefSMatthew G. Knepley   ierr = VecGetSize(X, &N);CHKERRQ(ierr);
1547*8e3a2eefSMatthew G. Knepley   ierr = MatSetSizes(*J, n, n, N, N);CHKERRQ(ierr);
1548*8e3a2eefSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1549*8e3a2eefSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) X);CHKERRQ(ierr);
1550*8e3a2eefSMatthew G. Knepley   ierr = PetscMalloc1(1, &ctx);CHKERRQ(ierr);
1551*8e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
1552*8e3a2eefSMatthew G. Knepley   ctx->X   = X;
1553*8e3a2eefSMatthew G. Knepley   ctx->ctx = user;
1554*8e3a2eefSMatthew G. Knepley   ierr = MatShellSetContext(*J, ctx);CHKERRQ(ierr);
1555*8e3a2eefSMatthew G. Knepley   ierr = MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void)) DMSNESJacobianMF_Destroy_Private);CHKERRQ(ierr);
1556*8e3a2eefSMatthew G. Knepley   ierr = MatShellSetOperation(*J, MATOP_MULT,    (void (*)(void)) DMSNESJacobianMF_Mult_Private);CHKERRQ(ierr);
1557*8e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
1558*8e3a2eefSMatthew G. Knepley }
1559*8e3a2eefSMatthew G. Knepley 
156038cfc46eSPierre Jolivet /*
156138cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
156238cfc46eSPierre Jolivet 
156338cfc46eSPierre Jolivet    Input Parameters:
156438cfc46eSPierre Jolivet +     X - SNES linearization point
156538cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
156638cfc46eSPierre Jolivet 
156738cfc46eSPierre Jolivet    Output Parameter:
156838cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
156938cfc46eSPierre Jolivet 
157038cfc46eSPierre Jolivet    Level: intermediate
157138cfc46eSPierre Jolivet 
157238cfc46eSPierre Jolivet .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
157338cfc46eSPierre Jolivet */
157438cfc46eSPierre Jolivet static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
157538cfc46eSPierre Jolivet {
157638cfc46eSPierre Jolivet   SNES           snes;
157738cfc46eSPierre Jolivet   Mat            pJ;
157838cfc46eSPierre Jolivet   DM             ovldm,origdm;
157938cfc46eSPierre Jolivet   DMSNES         sdm;
158038cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM,Vec,void*);
158138cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
158238cfc46eSPierre Jolivet   void           *bctx,*jctx;
158338cfc46eSPierre Jolivet   PetscErrorCode ierr;
158438cfc46eSPierre Jolivet 
158538cfc46eSPierre Jolivet   PetscFunctionBegin;
158638cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr);
158738cfc46eSPierre Jolivet   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr);
158838cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr);
158938cfc46eSPierre Jolivet   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr);
159038cfc46eSPierre Jolivet   ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr);
159138cfc46eSPierre Jolivet   ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr);
159238cfc46eSPierre Jolivet   ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr);
159338cfc46eSPierre Jolivet   ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr);
159438cfc46eSPierre Jolivet   ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr);
159538cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr);
159638cfc46eSPierre Jolivet   if (!snes) {
159738cfc46eSPierre Jolivet     ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr);
159838cfc46eSPierre Jolivet     ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr);
159938cfc46eSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr);
160038cfc46eSPierre Jolivet     ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr);
160138cfc46eSPierre Jolivet   }
160238cfc46eSPierre Jolivet   ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr);
160338cfc46eSPierre Jolivet   ierr = VecLockReadPush(X);CHKERRQ(ierr);
160438cfc46eSPierre Jolivet   PetscStackPush("SNES user Jacobian function");
160538cfc46eSPierre Jolivet   ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr);
160638cfc46eSPierre Jolivet   PetscStackPop;
160738cfc46eSPierre Jolivet   ierr = VecLockReadPop(X);CHKERRQ(ierr);
160838cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
160938cfc46eSPierre Jolivet   {
161038cfc46eSPierre Jolivet     Mat locpJ;
161138cfc46eSPierre Jolivet 
161238cfc46eSPierre Jolivet     ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr);
161338cfc46eSPierre Jolivet     ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
161438cfc46eSPierre Jolivet   }
161538cfc46eSPierre Jolivet   PetscFunctionReturn(0);
161638cfc46eSPierre Jolivet }
161738cfc46eSPierre Jolivet 
1618a925c78cSMatthew G. Knepley /*@
16199f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
16209f520fc2SToby Isaac 
16219f520fc2SToby Isaac   Input Parameters:
16229f520fc2SToby Isaac + dm - The DM object
1623dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
16249f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
16259f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
16261a244344SSatish Balay 
16271a244344SSatish Balay   Level: developer
16289f520fc2SToby Isaac @*/
16299f520fc2SToby Isaac PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
16309f520fc2SToby Isaac {
16319f520fc2SToby Isaac   PetscErrorCode ierr;
16329f520fc2SToby Isaac 
16339f520fc2SToby Isaac   PetscFunctionBegin;
16349f520fc2SToby Isaac   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
16359f520fc2SToby Isaac   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
16369f520fc2SToby Isaac   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
163738cfc46eSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr);
16389f520fc2SToby Isaac   PetscFunctionReturn(0);
16399f520fc2SToby Isaac }
16409f520fc2SToby Isaac 
1641f017bcb6SMatthew G. Knepley /*@C
1642f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1643f017bcb6SMatthew G. Knepley 
1644f017bcb6SMatthew G. Knepley   Input Parameters:
1645f017bcb6SMatthew G. Knepley + snes - the SNES object
1646f017bcb6SMatthew G. Knepley . dm   - the DM
1647f2cacb80SMatthew G. Knepley . t    - the time
1648f017bcb6SMatthew G. Knepley . u    - a DM vector
1649f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1650f017bcb6SMatthew G. Knepley 
1651f3c94560SMatthew G. Knepley   Output Parameters:
1652f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1653f3c94560SMatthew G. Knepley 
16547f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
16557f96f943SMatthew G. Knepley 
1656f017bcb6SMatthew G. Knepley   Level: developer
1657f017bcb6SMatthew G. Knepley 
16587f96f943SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1659f017bcb6SMatthew G. Knepley @*/
1660f2cacb80SMatthew G. Knepley PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
16611878804aSMatthew G. Knepley {
1662f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1663f017bcb6SMatthew G. Knepley   void            **ectxs;
1664f3c94560SMatthew G. Knepley   PetscReal        *err;
16657f96f943SMatthew G. Knepley   MPI_Comm          comm;
16667f96f943SMatthew G. Knepley   PetscInt          Nf, f;
16671878804aSMatthew G. Knepley   PetscErrorCode    ierr;
16681878804aSMatthew G. Knepley 
16691878804aSMatthew G. Knepley   PetscFunctionBegin;
1670f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1671f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1672064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1673534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
16747f96f943SMatthew G. Knepley 
1675f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
16767f96f943SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
16777f96f943SMatthew G. Knepley 
1678f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1679f017bcb6SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1680083401c6SMatthew G. Knepley   ierr = PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr);
16817f96f943SMatthew G. Knepley   {
16827f96f943SMatthew G. Knepley     PetscInt Nds, s;
16837f96f943SMatthew G. Knepley 
1684083401c6SMatthew G. Knepley     ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1685083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
16867f96f943SMatthew G. Knepley       PetscDS         ds;
1687083401c6SMatthew G. Knepley       DMLabel         label;
1688083401c6SMatthew G. Knepley       IS              fieldIS;
16897f96f943SMatthew G. Knepley       const PetscInt *fields;
16907f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1691083401c6SMatthew G. Knepley 
1692083401c6SMatthew G. Knepley       ierr = DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
1693083401c6SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
1694083401c6SMatthew G. Knepley       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
1695083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1696083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
1697083401c6SMatthew G. Knepley         ierr = PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);CHKERRQ(ierr);
1698083401c6SMatthew G. Knepley       }
1699083401c6SMatthew G. Knepley       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
1700083401c6SMatthew G. Knepley     }
1701083401c6SMatthew G. Knepley   }
1702f017bcb6SMatthew G. Knepley   if (Nf > 1) {
1703f2cacb80SMatthew G. Knepley     ierr = DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);CHKERRQ(ierr);
1704f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1705f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1706f3c94560SMatthew 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);
17071878804aSMatthew G. Knepley       }
1708b878532fSMatthew G. Knepley     } else if (error) {
1709b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
17101878804aSMatthew G. Knepley     } else {
1711f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
1712f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1713f017bcb6SMatthew G. Knepley         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
1714f3c94560SMatthew G. Knepley         ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr);
17151878804aSMatthew G. Knepley       }
1716f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
1717f017bcb6SMatthew G. Knepley     }
1718f017bcb6SMatthew G. Knepley   } else {
1719f2cacb80SMatthew G. Knepley     ierr = DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);CHKERRQ(ierr);
1720f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1721f3c94560SMatthew G. Knepley       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1722b878532fSMatthew G. Knepley     } else if (error) {
1723b878532fSMatthew G. Knepley       error[0] = err[0];
1724f017bcb6SMatthew G. Knepley     } else {
1725f3c94560SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);CHKERRQ(ierr);
1726f017bcb6SMatthew G. Knepley     }
1727f017bcb6SMatthew G. Knepley   }
1728f3c94560SMatthew G. Knepley   ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr);
1729f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1730f017bcb6SMatthew G. Knepley }
1731f017bcb6SMatthew G. Knepley 
1732f017bcb6SMatthew G. Knepley /*@C
1733f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1734f017bcb6SMatthew G. Knepley 
1735f017bcb6SMatthew G. Knepley   Input Parameters:
1736f017bcb6SMatthew G. Knepley + snes - the SNES object
1737f017bcb6SMatthew G. Knepley . dm   - the DM
1738f017bcb6SMatthew G. Knepley . u    - a DM vector
1739f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1740f017bcb6SMatthew G. Knepley 
1741f3c94560SMatthew G. Knepley   Output Parameters:
1742f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1743f3c94560SMatthew G. Knepley 
1744f017bcb6SMatthew G. Knepley   Level: developer
1745f017bcb6SMatthew G. Knepley 
1746f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1747f017bcb6SMatthew G. Knepley @*/
1748f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1749f017bcb6SMatthew G. Knepley {
1750f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1751f017bcb6SMatthew G. Knepley   Vec            r;
1752f017bcb6SMatthew G. Knepley   PetscReal      res;
1753f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1754f017bcb6SMatthew G. Knepley 
1755f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1756f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1757f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1758f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1759534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
1760f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1761f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1762f017bcb6SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
17631878804aSMatthew G. Knepley   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
17641878804aSMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
1765f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
1766f017bcb6SMatthew G. Knepley     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1767b878532fSMatthew G. Knepley   } else if (residual) {
1768b878532fSMatthew G. Knepley     *residual = res;
1769f017bcb6SMatthew G. Knepley   } else {
1770f017bcb6SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
17711878804aSMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
17721878804aSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
1773685405a1SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
1774685405a1SBarry Smith     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
1775f017bcb6SMatthew G. Knepley   }
1776f017bcb6SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
1777f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1778f017bcb6SMatthew G. Knepley }
1779f017bcb6SMatthew G. Knepley 
1780f017bcb6SMatthew G. Knepley /*@C
1781f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1782f017bcb6SMatthew G. Knepley 
1783f017bcb6SMatthew G. Knepley   Input Parameters:
1784f017bcb6SMatthew G. Knepley + snes - the SNES object
1785f017bcb6SMatthew G. Knepley . dm   - the DM
1786f017bcb6SMatthew G. Knepley . u    - a DM vector
1787f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1788f017bcb6SMatthew G. Knepley 
1789f3c94560SMatthew G. Knepley   Output Parameters:
1790f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1791f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1792f3c94560SMatthew G. Knepley 
1793f017bcb6SMatthew G. Knepley   Level: developer
1794f017bcb6SMatthew G. Knepley 
1795f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1796f017bcb6SMatthew G. Knepley @*/
1797f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1798f017bcb6SMatthew G. Knepley {
1799f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1800f017bcb6SMatthew G. Knepley   PetscDS        ds;
1801f017bcb6SMatthew G. Knepley   Mat            J, M;
1802f017bcb6SMatthew G. Knepley   MatNullSpace   nullspace;
1803f3c94560SMatthew G. Knepley   PetscReal      slope, intercept;
18047f96f943SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1805f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1806f017bcb6SMatthew G. Knepley 
1807f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1808f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1809f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1810f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1811534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1812064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 6);
1813f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1814f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1815f017bcb6SMatthew G. Knepley   /* Create and view matrices */
1816f017bcb6SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
18177f96f943SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1818f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1819f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1820282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
1821282e7bb4SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
1822282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
1823f017bcb6SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
1824282e7bb4SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
1825282e7bb4SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
1826282e7bb4SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
1827282e7bb4SMatthew G. Knepley   } else {
1828282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
1829282e7bb4SMatthew G. Knepley   }
1830f017bcb6SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
1831282e7bb4SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
1832282e7bb4SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
1833f017bcb6SMatthew G. Knepley   /* Check nullspace */
1834f017bcb6SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
1835f017bcb6SMatthew G. Knepley   if (nullspace) {
18361878804aSMatthew G. Knepley     PetscBool isNull;
1837f017bcb6SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
1838f017bcb6SMatthew G. Knepley     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
18391878804aSMatthew G. Knepley   }
1840f017bcb6SMatthew G. Knepley   /* Taylor test */
1841f017bcb6SMatthew G. Knepley   {
1842f017bcb6SMatthew G. Knepley     PetscRandom rand;
1843f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1844f3c94560SMatthew G. Knepley     PetscReal   h;
1845f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1846f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1847f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1848f017bcb6SMatthew G. Knepley 
1849f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
1850f017bcb6SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
1851f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
1852f017bcb6SMatthew G. Knepley     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
1853f017bcb6SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
1854f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
1855f017bcb6SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
1856f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
1857f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1858f017bcb6SMatthew G. Knepley     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1859f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1860f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1861f017bcb6SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
1862f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
1863f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
1864f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1865f017bcb6SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
1866f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1867f017bcb6SMatthew G. Knepley       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
1868f017bcb6SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
1869f017bcb6SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
1870f017bcb6SMatthew G. Knepley 
1871f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1872f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1873f017bcb6SMatthew G. Knepley     }
1874f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
1875f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
1876f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
18771878804aSMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
1878f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
1879f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1880f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1881f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
1882f017bcb6SMatthew G. Knepley     }
1883f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
1884f017bcb6SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
1885f017bcb6SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
1886f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1887f017bcb6SMatthew G. Knepley     if (tol >= 0) {
1888b878532fSMatthew 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);
1889b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1890b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1891b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1892f017bcb6SMatthew G. Knepley     } else {
1893f3c94560SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
1894f017bcb6SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
1895f017bcb6SMatthew G. Knepley     }
1896f017bcb6SMatthew G. Knepley   }
18971878804aSMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
18981878804aSMatthew G. Knepley   PetscFunctionReturn(0);
18991878804aSMatthew G. Knepley }
19001878804aSMatthew G. Knepley 
19017f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1902f017bcb6SMatthew G. Knepley {
1903f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1904f017bcb6SMatthew G. Knepley 
1905f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1906f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);CHKERRQ(ierr);
1907f3c94560SMatthew G. Knepley   ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr);
1908f3c94560SMatthew G. Knepley   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr);
1909f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1910f017bcb6SMatthew G. Knepley }
1911f017bcb6SMatthew G. Knepley 
1912bee9a294SMatthew G. Knepley /*@C
1913bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1914bee9a294SMatthew G. Knepley 
1915bee9a294SMatthew G. Knepley   Input Parameters:
1916bee9a294SMatthew G. Knepley + snes - the SNES object
19177f96f943SMatthew G. Knepley - u    - representative SNES vector
19187f96f943SMatthew G. Knepley 
19197f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1920bee9a294SMatthew G. Knepley 
1921bee9a294SMatthew G. Knepley   Level: developer
1922bee9a294SMatthew G. Knepley @*/
19237f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
19241878804aSMatthew G. Knepley {
19251878804aSMatthew G. Knepley   DM             dm;
19261878804aSMatthew G. Knepley   Vec            sol;
19271878804aSMatthew G. Knepley   PetscBool      check;
19281878804aSMatthew G. Knepley   PetscErrorCode ierr;
19291878804aSMatthew G. Knepley 
19301878804aSMatthew G. Knepley   PetscFunctionBegin;
1931c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
19321878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
19331878804aSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
19341878804aSMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
19351878804aSMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
19367f96f943SMatthew G. Knepley   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
19371878804aSMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1938552f7358SJed Brown   PetscFunctionReturn(0);
1939552f7358SJed Brown }
1940