xref: /petsc/src/snes/utils/dmplexsnes.c (revision 680d18d566775cd65817f8ecb7f0a4c7614c2523)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>     /*I "petscsnes.h"   I*/
324cdb843SMatthew G. Knepley #include <petscds.h>
4a925c78cSMatthew G. Knepley #include <petscblaslapack.h>
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
6af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
7552f7358SJed Brown 
8649ef022SMatthew Knepley static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
9649ef022SMatthew Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
10649ef022SMatthew Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
11649ef022SMatthew Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
12649ef022SMatthew Knepley {
13649ef022SMatthew Knepley   p[0] = u[uOff[1]];
14649ef022SMatthew Knepley }
15649ef022SMatthew Knepley 
16649ef022SMatthew Knepley /*
17649ef022SMatthew Knepley   SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately.
18649ef022SMatthew Knepley 
19649ef022SMatthew Knepley   Collective on SNES
20649ef022SMatthew Knepley 
21649ef022SMatthew Knepley   Input Parameters:
22649ef022SMatthew Knepley + snes      - The SNES
23649ef022SMatthew Knepley . pfield    - The field number for pressure
24649ef022SMatthew Knepley . nullspace - The pressure nullspace
25649ef022SMatthew Knepley . u         - The solution vector
26649ef022SMatthew Knepley - ctx       - An optional user context
27649ef022SMatthew Knepley 
28649ef022SMatthew Knepley   Output Parameter:
29649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
30649ef022SMatthew Knepley 
31649ef022SMatthew Knepley   Notes:
32649ef022SMatthew Knepley   If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly.
33649ef022SMatthew Knepley 
34649ef022SMatthew Knepley   Level: developer
35649ef022SMatthew Knepley 
36649ef022SMatthew Knepley .seealso: SNESConvergedCorrectPressure()
37649ef022SMatthew Knepley */
38649ef022SMatthew Knepley static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
39649ef022SMatthew Knepley {
40649ef022SMatthew Knepley   DM             dm;
41649ef022SMatthew Knepley   PetscDS        ds;
42649ef022SMatthew Knepley   const Vec     *nullvecs;
43649ef022SMatthew Knepley   PetscScalar    pintd, *intc, *intn;
44649ef022SMatthew Knepley   MPI_Comm       comm;
45649ef022SMatthew Knepley   PetscInt       Nf, Nv;
46649ef022SMatthew Knepley   PetscErrorCode ierr;
47649ef022SMatthew Knepley 
48649ef022SMatthew Knepley   PetscFunctionBegin;
49649ef022SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
50649ef022SMatthew Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
51649ef022SMatthew Knepley   if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
52649ef022SMatthew Knepley   if (!nullspace) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
53649ef022SMatthew Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
54649ef022SMatthew Knepley   ierr = PetscDSSetObjective(ds, pfield, pressure_Private);CHKERRQ(ierr);
55649ef022SMatthew Knepley   ierr = MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs);CHKERRQ(ierr);
56649ef022SMatthew Knepley   if (Nv != 1) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %D", Nv);
57649ef022SMatthew Knepley   ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
58649ef022SMatthew Knepley   if (PetscAbsScalar(pintd) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
59649ef022SMatthew Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
60649ef022SMatthew Knepley   ierr = PetscMalloc2(Nf, &intc, Nf, &intn);CHKERRQ(ierr);
61649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx);CHKERRQ(ierr);
62649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
63649ef022SMatthew Knepley   ierr = VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]);CHKERRQ(ierr);
64649ef022SMatthew Knepley #if defined (PETSC_USE_DEBUG)
65649ef022SMatthew Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
66649ef022SMatthew Knepley   if (PetscAbsScalar(intc[pfield]) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g\n", (double) PetscRealPart(intc[pfield]));
67649ef022SMatthew Knepley #endif
68649ef022SMatthew Knepley   ierr = PetscFree2(intc, intn);CHKERRQ(ierr);
69649ef022SMatthew Knepley   PetscFunctionReturn(0);
70649ef022SMatthew Knepley }
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley /*@C
73649ef022SMatthew Knepley    SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics SNESConvergedDefault().
74649ef022SMatthew Knepley 
75649ef022SMatthew Knepley    Logically Collective on SNES
76649ef022SMatthew Knepley 
77649ef022SMatthew Knepley    Input Parameters:
78649ef022SMatthew Knepley +  snes - the SNES context
79649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
80649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
81649ef022SMatthew Knepley .  snorm - 2-norm of current step
82649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
83649ef022SMatthew Knepley -  ctx   - Optional user context
84649ef022SMatthew Knepley 
85649ef022SMatthew Knepley    Output Parameter:
86649ef022SMatthew Knepley .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
87649ef022SMatthew Knepley 
88649ef022SMatthew Knepley    Notes:
89649ef022SMatthew Knepley    In order to use this monitor, you must setup several PETSc structures. First fields must be added to the DM, and a PetscDS must be created with discretizations of those fields. We currently assume that the pressure field has index 1. The pressure field must have a nullspace, likely created using the DMSetNullSpaceConstructor() interface. Last we must be able to integrate the pressure over the domain, so the DM attached to the SNES must be a Plex at this time.
90649ef022SMatthew Knepley 
91649ef022SMatthew Knepley    Level: advanced
92649ef022SMatthew Knepley 
93649ef022SMatthew Knepley .seealso: SNESConvergedDefault(), SNESSetConvergenceTest(), DMSetNullSpaceConstructor()
94649ef022SMatthew Knepley @*/
95649ef022SMatthew Knepley PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
96649ef022SMatthew Knepley {
97649ef022SMatthew Knepley   PetscBool      monitorIntegral = PETSC_FALSE;
98649ef022SMatthew Knepley   PetscErrorCode ierr;
99649ef022SMatthew Knepley 
100649ef022SMatthew Knepley   PetscFunctionBegin;
101649ef022SMatthew Knepley   ierr = SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx);CHKERRQ(ierr);
102649ef022SMatthew Knepley   if (monitorIntegral) {
103649ef022SMatthew Knepley     Mat          J;
104649ef022SMatthew Knepley     Vec          u;
105649ef022SMatthew Knepley     MatNullSpace nullspace;
106649ef022SMatthew Knepley     const Vec   *nullvecs;
107649ef022SMatthew Knepley     PetscScalar  pintd;
108649ef022SMatthew Knepley 
109649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
110649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
111649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
112649ef022SMatthew Knepley     ierr = MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);CHKERRQ(ierr);
113649ef022SMatthew Knepley     ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
114649ef022SMatthew Knepley     ierr = PetscInfo1(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));CHKERRQ(ierr);
115649ef022SMatthew Knepley   }
116649ef022SMatthew Knepley   if (*reason > 0) {
117649ef022SMatthew Knepley     Mat          J;
118649ef022SMatthew Knepley     Vec          u;
119649ef022SMatthew Knepley     MatNullSpace nullspace;
120649ef022SMatthew Knepley     PetscInt     pfield = 1;
121649ef022SMatthew Knepley 
122649ef022SMatthew Knepley     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
123649ef022SMatthew Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
124649ef022SMatthew Knepley     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
125649ef022SMatthew Knepley     ierr = SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx);CHKERRQ(ierr);
126649ef022SMatthew Knepley   }
127649ef022SMatthew Knepley   PetscFunctionReturn(0);
128649ef022SMatthew Knepley }
129649ef022SMatthew Knepley 
13024cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
13124cdb843SMatthew G. Knepley 
1326da023fcSToby Isaac static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
1336da023fcSToby Isaac {
1346da023fcSToby Isaac   PetscBool      isPlex;
1356da023fcSToby Isaac   PetscErrorCode ierr;
1366da023fcSToby Isaac 
1376da023fcSToby Isaac   PetscFunctionBegin;
1386da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
1396da023fcSToby Isaac   if (isPlex) {
1406da023fcSToby Isaac     *plex = dm;
1416da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
142f7148743SMatthew G. Knepley   } else {
143f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
144f7148743SMatthew G. Knepley     if (!*plex) {
1456da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
146f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
1476da023fcSToby Isaac       if (copy) {
1486da023fcSToby Isaac         PetscInt    i;
1496da023fcSToby Isaac         PetscObject obj;
1506da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
1516da023fcSToby Isaac 
1526da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
1536da023fcSToby Isaac         for (i = 0; i < 3; i++) {
1546da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
1556da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
1566da023fcSToby Isaac         }
1576da023fcSToby Isaac       }
158f7148743SMatthew G. Knepley     } else {
159f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
160f7148743SMatthew G. Knepley     }
1616da023fcSToby Isaac   }
1626da023fcSToby Isaac   PetscFunctionReturn(0);
1636da023fcSToby Isaac }
1646da023fcSToby Isaac 
1654267b1a3SMatthew G. Knepley /*@C
1664267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
1674267b1a3SMatthew G. Knepley 
168d083f849SBarry Smith   Collective
1694267b1a3SMatthew G. Knepley 
1704267b1a3SMatthew G. Knepley   Input Parameter:
1714267b1a3SMatthew G. Knepley . comm - the communicator
1724267b1a3SMatthew G. Knepley 
1734267b1a3SMatthew G. Knepley   Output Parameter:
1744267b1a3SMatthew G. Knepley . ctx - the context
1754267b1a3SMatthew G. Knepley 
1764267b1a3SMatthew G. Knepley   Level: beginner
1774267b1a3SMatthew G. Knepley 
1784267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
1794267b1a3SMatthew G. Knepley @*/
1800adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
1810adebc6cSBarry Smith {
182552f7358SJed Brown   PetscErrorCode ierr;
183552f7358SJed Brown 
184552f7358SJed Brown   PetscFunctionBegin;
185552f7358SJed Brown   PetscValidPointer(ctx, 2);
18695dccacaSBarry Smith   ierr = PetscNew(ctx);CHKERRQ(ierr);
1871aa26658SKarl Rupp 
188552f7358SJed Brown   (*ctx)->comm   = comm;
189552f7358SJed Brown   (*ctx)->dim    = -1;
190552f7358SJed Brown   (*ctx)->nInput = 0;
1910298fd71SBarry Smith   (*ctx)->points = NULL;
1920298fd71SBarry Smith   (*ctx)->cells  = NULL;
193552f7358SJed Brown   (*ctx)->n      = -1;
1940298fd71SBarry Smith   (*ctx)->coords = NULL;
195552f7358SJed Brown   PetscFunctionReturn(0);
196552f7358SJed Brown }
197552f7358SJed Brown 
1984267b1a3SMatthew G. Knepley /*@C
1994267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
2004267b1a3SMatthew G. Knepley 
2014267b1a3SMatthew G. Knepley   Not collective
2024267b1a3SMatthew G. Knepley 
2034267b1a3SMatthew G. Knepley   Input Parameters:
2044267b1a3SMatthew G. Knepley + ctx - the context
2054267b1a3SMatthew G. Knepley - dim - the spatial dimension
2064267b1a3SMatthew G. Knepley 
2074267b1a3SMatthew G. Knepley   Level: intermediate
2084267b1a3SMatthew G. Knepley 
2094267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2104267b1a3SMatthew G. Knepley @*/
2110adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
2120adebc6cSBarry Smith {
213552f7358SJed Brown   PetscFunctionBegin;
214ff1e0c32SBarry Smith   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
215552f7358SJed Brown   ctx->dim = dim;
216552f7358SJed Brown   PetscFunctionReturn(0);
217552f7358SJed Brown }
218552f7358SJed Brown 
2194267b1a3SMatthew G. Knepley /*@C
2204267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2214267b1a3SMatthew G. Knepley 
2224267b1a3SMatthew G. Knepley   Not collective
2234267b1a3SMatthew G. Knepley 
2244267b1a3SMatthew G. Knepley   Input Parameter:
2254267b1a3SMatthew G. Knepley . ctx - the context
2264267b1a3SMatthew G. Knepley 
2274267b1a3SMatthew G. Knepley   Output Parameter:
2284267b1a3SMatthew G. Knepley . dim - the spatial dimension
2294267b1a3SMatthew G. Knepley 
2304267b1a3SMatthew G. Knepley   Level: intermediate
2314267b1a3SMatthew G. Knepley 
2324267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2334267b1a3SMatthew G. Knepley @*/
2340adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
2350adebc6cSBarry Smith {
236552f7358SJed Brown   PetscFunctionBegin;
237552f7358SJed Brown   PetscValidIntPointer(dim, 2);
238552f7358SJed Brown   *dim = ctx->dim;
239552f7358SJed Brown   PetscFunctionReturn(0);
240552f7358SJed Brown }
241552f7358SJed Brown 
2424267b1a3SMatthew G. Knepley /*@C
2434267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2444267b1a3SMatthew G. Knepley 
2454267b1a3SMatthew G. Knepley   Not collective
2464267b1a3SMatthew G. Knepley 
2474267b1a3SMatthew G. Knepley   Input Parameters:
2484267b1a3SMatthew G. Knepley + ctx - the context
2494267b1a3SMatthew G. Knepley - dof - the number of fields
2504267b1a3SMatthew G. Knepley 
2514267b1a3SMatthew G. Knepley   Level: intermediate
2524267b1a3SMatthew G. Knepley 
2534267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2544267b1a3SMatthew G. Knepley @*/
2550adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
2560adebc6cSBarry Smith {
257552f7358SJed Brown   PetscFunctionBegin;
258ff1e0c32SBarry Smith   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
259552f7358SJed Brown   ctx->dof = dof;
260552f7358SJed Brown   PetscFunctionReturn(0);
261552f7358SJed Brown }
262552f7358SJed Brown 
2634267b1a3SMatthew G. Knepley /*@C
2644267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2654267b1a3SMatthew G. Knepley 
2664267b1a3SMatthew G. Knepley   Not collective
2674267b1a3SMatthew G. Knepley 
2684267b1a3SMatthew G. Knepley   Input Parameter:
2694267b1a3SMatthew G. Knepley . ctx - the context
2704267b1a3SMatthew G. Knepley 
2714267b1a3SMatthew G. Knepley   Output Parameter:
2724267b1a3SMatthew G. Knepley . dof - the number of fields
2734267b1a3SMatthew G. Knepley 
2744267b1a3SMatthew G. Knepley   Level: intermediate
2754267b1a3SMatthew G. Knepley 
2764267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
2774267b1a3SMatthew G. Knepley @*/
2780adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
2790adebc6cSBarry Smith {
280552f7358SJed Brown   PetscFunctionBegin;
281552f7358SJed Brown   PetscValidIntPointer(dof, 2);
282552f7358SJed Brown   *dof = ctx->dof;
283552f7358SJed Brown   PetscFunctionReturn(0);
284552f7358SJed Brown }
285552f7358SJed Brown 
2864267b1a3SMatthew G. Knepley /*@C
2874267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2884267b1a3SMatthew G. Knepley 
2894267b1a3SMatthew G. Knepley   Not collective
2904267b1a3SMatthew G. Knepley 
2914267b1a3SMatthew G. Knepley   Input Parameters:
2924267b1a3SMatthew G. Knepley + ctx    - the context
2934267b1a3SMatthew G. Knepley . n      - the number of points
2944267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2954267b1a3SMatthew G. Knepley 
2964267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
2974267b1a3SMatthew G. Knepley 
2984267b1a3SMatthew G. Knepley   Level: intermediate
2994267b1a3SMatthew G. Knepley 
3004267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
3014267b1a3SMatthew G. Knepley @*/
3020adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
3030adebc6cSBarry Smith {
304552f7358SJed Brown   PetscErrorCode ierr;
305552f7358SJed Brown 
306552f7358SJed Brown   PetscFunctionBegin;
3070adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
3080adebc6cSBarry Smith   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
309552f7358SJed Brown   ctx->nInput = n;
3101aa26658SKarl Rupp 
311785e854fSJed Brown   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
312580bdb30SBarry Smith   ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr);
313552f7358SJed Brown   PetscFunctionReturn(0);
314552f7358SJed Brown }
315552f7358SJed Brown 
3164267b1a3SMatthew G. Knepley /*@C
3174267b1a3SMatthew G. Knepley   DMInterpolationSetUp - Computea spatial indices that add in point location during interpolation
3184267b1a3SMatthew G. Knepley 
3194267b1a3SMatthew G. Knepley   Collective on ctx
3204267b1a3SMatthew G. Knepley 
3214267b1a3SMatthew G. Knepley   Input Parameters:
3224267b1a3SMatthew G. Knepley + ctx - the context
3234267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
3244267b1a3SMatthew G. Knepley - redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
3254267b1a3SMatthew G. Knepley 
3264267b1a3SMatthew G. Knepley   Level: intermediate
3274267b1a3SMatthew G. Knepley 
3284267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3294267b1a3SMatthew G. Knepley @*/
3300adebc6cSBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
3310adebc6cSBarry Smith {
332552f7358SJed Brown   MPI_Comm          comm = ctx->comm;
333552f7358SJed Brown   PetscScalar       *a;
334552f7358SJed Brown   PetscInt          p, q, i;
335552f7358SJed Brown   PetscMPIInt       rank, size;
336552f7358SJed Brown   PetscErrorCode    ierr;
337552f7358SJed Brown   Vec               pointVec;
3383a93e3b7SToby Isaac   PetscSF           cellSF;
339552f7358SJed Brown   PetscLayout       layout;
340552f7358SJed Brown   PetscReal         *globalPoints;
341cb313848SJed Brown   PetscScalar       *globalPointsScalar;
342552f7358SJed Brown   const PetscInt    *ranges;
343552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3443a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3453a93e3b7SToby Isaac   const PetscInt    *foundPoints;
346552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3473a93e3b7SToby Isaac   PetscInt          n, N, numFound;
348552f7358SJed Brown 
34919436ca2SJed Brown   PetscFunctionBegin;
35019436ca2SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
351ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
352ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
3530adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
35419436ca2SJed Brown   /* Locate points */
35519436ca2SJed Brown   n = ctx->nInput;
356552f7358SJed Brown   if (!redundantPoints) {
357552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
358552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
359552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
360552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
361552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
362552f7358SJed Brown     /* Communicate all points to all processes */
363dcca6d9dSJed Brown     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
364552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
365552f7358SJed Brown     for (p = 0; p < size; ++p) {
366552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
367552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
368552f7358SJed Brown     }
369ffc4695bSBarry Smith     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRMPI(ierr);
370552f7358SJed Brown   } else {
371552f7358SJed Brown     N = n;
372552f7358SJed Brown     globalPoints = ctx->points;
37338ea73c8SJed Brown     counts = displs = NULL;
37438ea73c8SJed Brown     layout = NULL;
375552f7358SJed Brown   }
376552f7358SJed Brown #if 0
377dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
37819436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
379552f7358SJed Brown #else
380cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
38146b3086cSToby Isaac   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
38246b3086cSToby Isaac   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
383cb313848SJed Brown #else
384cb313848SJed Brown   globalPointsScalar = globalPoints;
385cb313848SJed Brown #endif
38604706141SMatthew G Knepley   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
387dcca6d9dSJed Brown   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
3887b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
3893a93e3b7SToby Isaac   cellSF = NULL;
3902d1fa6caSMatthew G. Knepley   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
3913a93e3b7SToby Isaac   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
392552f7358SJed Brown #endif
3933a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3943a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
395552f7358SJed Brown   }
396552f7358SJed Brown   /* Let the lowest rank process own each point */
397b2566f29SBarry Smith   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
398552f7358SJed Brown   ctx->n = 0;
399552f7358SJed Brown   for (p = 0; p < N; ++p) {
400e5b268a4SToby Isaac     if (globalProcs[p] == size) 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));
4011aa26658SKarl Rupp     else if (globalProcs[p] == rank) ctx->n++;
402552f7358SJed Brown   }
403552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
404785e854fSJed Brown   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
405552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
406552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
407552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
408c0dedaeaSBarry Smith   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
409552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
410552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
411552f7358SJed Brown     if (globalProcs[p] == rank) {
412552f7358SJed Brown       PetscInt d;
413552f7358SJed Brown 
4141aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
415f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
416f317b747SMatthew G. Knepley       ++q;
417552f7358SJed Brown     }
418552f7358SJed Brown   }
419552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
420552f7358SJed Brown #if 0
421552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
422552f7358SJed Brown #else
423552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
4243a93e3b7SToby Isaac   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
425552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
426552f7358SJed Brown #endif
427cb313848SJed Brown   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
428d343d804SMatthew G. Knepley   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
429552f7358SJed Brown   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
430552f7358SJed Brown   PetscFunctionReturn(0);
431552f7358SJed Brown }
432552f7358SJed Brown 
4334267b1a3SMatthew G. Knepley /*@C
4344267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
4354267b1a3SMatthew G. Knepley 
4364267b1a3SMatthew G. Knepley   Collective on ctx
4374267b1a3SMatthew G. Knepley 
4384267b1a3SMatthew G. Knepley   Input Parameter:
4394267b1a3SMatthew G. Knepley . ctx - the context
4404267b1a3SMatthew G. Knepley 
4414267b1a3SMatthew G. Knepley   Output Parameter:
4424267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4434267b1a3SMatthew G. Knepley 
4444267b1a3SMatthew 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.
4454267b1a3SMatthew G. Knepley 
4464267b1a3SMatthew G. Knepley   Level: intermediate
4474267b1a3SMatthew G. Knepley 
4484267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4494267b1a3SMatthew G. Knepley @*/
4500adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
4510adebc6cSBarry Smith {
452552f7358SJed Brown   PetscFunctionBegin;
453552f7358SJed Brown   PetscValidPointer(coordinates, 2);
4540adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
455552f7358SJed Brown   *coordinates = ctx->coords;
456552f7358SJed Brown   PetscFunctionReturn(0);
457552f7358SJed Brown }
458552f7358SJed Brown 
4594267b1a3SMatthew G. Knepley /*@C
4604267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
4614267b1a3SMatthew G. Knepley 
4624267b1a3SMatthew G. Knepley   Collective on ctx
4634267b1a3SMatthew G. Knepley 
4644267b1a3SMatthew G. Knepley   Input Parameter:
4654267b1a3SMatthew G. Knepley . ctx - the context
4664267b1a3SMatthew G. Knepley 
4674267b1a3SMatthew G. Knepley   Output Parameter:
4684267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4694267b1a3SMatthew G. Knepley 
4704267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
4714267b1a3SMatthew G. Knepley 
4724267b1a3SMatthew G. Knepley   Level: intermediate
4734267b1a3SMatthew G. Knepley 
4744267b1a3SMatthew G. Knepley .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
4754267b1a3SMatthew G. Knepley @*/
4760adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
4770adebc6cSBarry Smith {
478552f7358SJed Brown   PetscErrorCode ierr;
479552f7358SJed Brown 
480552f7358SJed Brown   PetscFunctionBegin;
481552f7358SJed Brown   PetscValidPointer(v, 2);
4820adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
483552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
484552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
485552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
486c0dedaeaSBarry Smith   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
487552f7358SJed Brown   PetscFunctionReturn(0);
488552f7358SJed Brown }
489552f7358SJed Brown 
4904267b1a3SMatthew G. Knepley /*@C
4914267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
4924267b1a3SMatthew G. Knepley 
4934267b1a3SMatthew G. Knepley   Collective on ctx
4944267b1a3SMatthew G. Knepley 
4954267b1a3SMatthew G. Knepley   Input Parameters:
4964267b1a3SMatthew G. Knepley + ctx - the context
4974267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
4984267b1a3SMatthew G. Knepley 
4994267b1a3SMatthew G. Knepley   Level: intermediate
5004267b1a3SMatthew G. Knepley 
5014267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
5024267b1a3SMatthew G. Knepley @*/
5030adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
5040adebc6cSBarry Smith {
505552f7358SJed Brown   PetscErrorCode ierr;
506552f7358SJed Brown 
507552f7358SJed Brown   PetscFunctionBegin;
508552f7358SJed Brown   PetscValidPointer(v, 2);
5090adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
510552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
511552f7358SJed Brown   PetscFunctionReturn(0);
512552f7358SJed Brown }
513552f7358SJed Brown 
5147a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
515a6dfd86eSKarl Rupp {
516552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
51756044e6dSMatthew G. Knepley   const PetscScalar *coords;
51856044e6dSMatthew G. Knepley   PetscScalar    *a;
519552f7358SJed Brown   PetscInt       p;
520552f7358SJed Brown   PetscErrorCode ierr;
521552f7358SJed Brown 
522552f7358SJed Brown   PetscFunctionBegin;
523dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
52456044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
525552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
526552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
527552f7358SJed Brown     PetscInt     c = ctx->cells[p];
528a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
529552f7358SJed Brown     PetscReal    xi[4];
530552f7358SJed Brown     PetscInt     d, f, comp;
531552f7358SJed Brown 
5328e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
533ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5340298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5351aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5361aa26658SKarl Rupp 
537552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
538552f7358SJed Brown       xi[d] = 0.0;
5391aa26658SKarl 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]);
5401aa26658SKarl 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];
541552f7358SJed Brown     }
5420298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
543552f7358SJed Brown   }
544552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
54556044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
546552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
547552f7358SJed Brown   PetscFunctionReturn(0);
548552f7358SJed Brown }
549552f7358SJed Brown 
5507a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
5517a1931ceSMatthew G. Knepley {
5527a1931ceSMatthew G. Knepley   PetscReal      *v0, *J, *invJ, detJ;
55356044e6dSMatthew G. Knepley   const PetscScalar *coords;
55456044e6dSMatthew G. Knepley   PetscScalar    *a;
5557a1931ceSMatthew G. Knepley   PetscInt       p;
5567a1931ceSMatthew G. Knepley   PetscErrorCode ierr;
5577a1931ceSMatthew G. Knepley 
5587a1931ceSMatthew G. Knepley   PetscFunctionBegin;
559dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
56056044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5617a1931ceSMatthew G. Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
5627a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5637a1931ceSMatthew G. Knepley     PetscInt       c = ctx->cells[p];
5647a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5652584bbe8SMatthew G. Knepley     PetscScalar   *x = NULL;
5667a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5677a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5687a1931ceSMatthew G. Knepley 
5698e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
570ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
5717a1931ceSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5727a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5737a1931ceSMatthew G. Knepley 
5747a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5757a1931ceSMatthew G. Knepley       xi[d] = 0.0;
5767a1931ceSMatthew 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]);
5777a1931ceSMatthew 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];
5787a1931ceSMatthew G. Knepley     }
5797a1931ceSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
5807a1931ceSMatthew G. Knepley   }
5817a1931ceSMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
58256044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
5837a1931ceSMatthew G. Knepley   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
5847a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
5857a1931ceSMatthew G. Knepley }
5867a1931ceSMatthew G. Knepley 
5875820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
588552f7358SJed Brown {
589552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
590552f7358SJed Brown   const PetscScalar x0        = vertices[0];
591552f7358SJed Brown   const PetscScalar y0        = vertices[1];
592552f7358SJed Brown   const PetscScalar x1        = vertices[2];
593552f7358SJed Brown   const PetscScalar y1        = vertices[3];
594552f7358SJed Brown   const PetscScalar x2        = vertices[4];
595552f7358SJed Brown   const PetscScalar y2        = vertices[5];
596552f7358SJed Brown   const PetscScalar x3        = vertices[6];
597552f7358SJed Brown   const PetscScalar y3        = vertices[7];
598552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
599552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
600552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
601552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
602552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
603552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
60456044e6dSMatthew G. Knepley   const PetscScalar *ref;
60556044e6dSMatthew G. Knepley   PetscScalar       *real;
606552f7358SJed Brown   PetscErrorCode    ierr;
607552f7358SJed Brown 
608552f7358SJed Brown   PetscFunctionBegin;
60956044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
610552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
611552f7358SJed Brown   {
612552f7358SJed Brown     const PetscScalar p0 = ref[0];
613552f7358SJed Brown     const PetscScalar p1 = ref[1];
614552f7358SJed Brown 
615552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
616552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
617552f7358SJed Brown   }
618552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
61956044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
620552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
621552f7358SJed Brown   PetscFunctionReturn(0);
622552f7358SJed Brown }
623552f7358SJed Brown 
624af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
625d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
626552f7358SJed Brown {
627552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
628552f7358SJed Brown   const PetscScalar x0        = vertices[0];
629552f7358SJed Brown   const PetscScalar y0        = vertices[1];
630552f7358SJed Brown   const PetscScalar x1        = vertices[2];
631552f7358SJed Brown   const PetscScalar y1        = vertices[3];
632552f7358SJed Brown   const PetscScalar x2        = vertices[4];
633552f7358SJed Brown   const PetscScalar y2        = vertices[5];
634552f7358SJed Brown   const PetscScalar x3        = vertices[6];
635552f7358SJed Brown   const PetscScalar y3        = vertices[7];
636552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
637552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
63856044e6dSMatthew G. Knepley   const PetscScalar *ref;
639552f7358SJed Brown   PetscErrorCode    ierr;
640552f7358SJed Brown 
641552f7358SJed Brown   PetscFunctionBegin;
64256044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
643552f7358SJed Brown   {
644552f7358SJed Brown     const PetscScalar x       = ref[0];
645552f7358SJed Brown     const PetscScalar y       = ref[1];
646552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
647da80777bSKarl Rupp     PetscScalar       values[4];
648da80777bSKarl Rupp 
649da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
650da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
65194ab13aaSBarry Smith     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
652552f7358SJed Brown   }
653552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
65456044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
65594ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65694ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657552f7358SJed Brown   PetscFunctionReturn(0);
658552f7358SJed Brown }
659552f7358SJed Brown 
660a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
661a6dfd86eSKarl Rupp {
662fafc0619SMatthew G Knepley   DM             dmCoord;
663de73a395SMatthew G. Knepley   PetscFE        fem = NULL;
664552f7358SJed Brown   SNES           snes;
665552f7358SJed Brown   KSP            ksp;
666552f7358SJed Brown   PC             pc;
667552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
668552f7358SJed Brown   Mat            J;
669ef0bb6c7SMatthew G. Knepley   PetscTabulation    T;
67056044e6dSMatthew G. Knepley   const PetscScalar *coords;
67156044e6dSMatthew G. Knepley   PetscScalar    *a;
672ef0bb6c7SMatthew G. Knepley   PetscReal      xir[2];
673de73a395SMatthew G. Knepley   PetscInt       Nf, p;
6745509d985SMatthew G. Knepley   const PetscInt dof = ctx->dof;
675552f7358SJed Brown   PetscErrorCode ierr;
676552f7358SJed Brown 
677552f7358SJed Brown   PetscFunctionBegin;
678de73a395SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
67944a7f3ddSMatthew G. Knepley   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
680552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
681fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
682552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
683552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
684552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
685552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
686c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
687552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
688552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
689552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
690552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
691552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
692552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
6930298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
6940298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
695552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
696552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
697552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
698552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
699552f7358SJed Brown 
70056044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
701552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
702ef0bb6c7SMatthew G. Knepley   ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr);
703552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
704a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
705552f7358SJed Brown     PetscScalar *xi;
706552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
707552f7358SJed Brown 
708552f7358SJed Brown     /* Can make this do all points at once */
7090298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
710ff1e0c32SBarry Smith     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
7110298fd71SBarry Smith     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
7120298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
7130298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
714552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
715552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
716552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
717552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
718552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
719552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
720cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
721cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
7225509d985SMatthew G. Knepley     if (4*dof != xSize) {
7235509d985SMatthew G. Knepley       PetscInt d;
7241aa26658SKarl Rupp 
7255509d985SMatthew G. Knepley       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
726ef0bb6c7SMatthew G. Knepley       ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr);
7275509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7285509d985SMatthew G. Knepley         a[p*dof+comp] = 0.0;
7295509d985SMatthew G. Knepley         for (d = 0; d < xSize/dof; ++d) {
730ef0bb6c7SMatthew G. Knepley           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
7315509d985SMatthew G. Knepley         }
7325509d985SMatthew G. Knepley       }
7335509d985SMatthew G. Knepley     } else {
7345509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp)
7355509d985SMatthew 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];
7365509d985SMatthew G. Knepley     }
737552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
7380298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
7390298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
740552f7358SJed Brown   }
741ef0bb6c7SMatthew G. Knepley   ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
742552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
74356044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
744552f7358SJed Brown 
745552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
746552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
747552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
748552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
749552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
750552f7358SJed Brown   PetscFunctionReturn(0);
751552f7358SJed Brown }
752552f7358SJed Brown 
7535820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
754552f7358SJed Brown {
755552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
756552f7358SJed Brown   const PetscScalar x0        = vertices[0];
757552f7358SJed Brown   const PetscScalar y0        = vertices[1];
758552f7358SJed Brown   const PetscScalar z0        = vertices[2];
7597a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
7607a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
7617a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
762552f7358SJed Brown   const PetscScalar x2        = vertices[6];
763552f7358SJed Brown   const PetscScalar y2        = vertices[7];
764552f7358SJed Brown   const PetscScalar z2        = vertices[8];
7657a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
7667a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
7677a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
768552f7358SJed Brown   const PetscScalar x4        = vertices[12];
769552f7358SJed Brown   const PetscScalar y4        = vertices[13];
770552f7358SJed Brown   const PetscScalar z4        = vertices[14];
771552f7358SJed Brown   const PetscScalar x5        = vertices[15];
772552f7358SJed Brown   const PetscScalar y5        = vertices[16];
773552f7358SJed Brown   const PetscScalar z5        = vertices[17];
774552f7358SJed Brown   const PetscScalar x6        = vertices[18];
775552f7358SJed Brown   const PetscScalar y6        = vertices[19];
776552f7358SJed Brown   const PetscScalar z6        = vertices[20];
777552f7358SJed Brown   const PetscScalar x7        = vertices[21];
778552f7358SJed Brown   const PetscScalar y7        = vertices[22];
779552f7358SJed Brown   const PetscScalar z7        = vertices[23];
780552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
781552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
782552f7358SJed Brown   const PetscScalar h_1       = z1 - z0;
783552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
784552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
785552f7358SJed Brown   const PetscScalar h_3       = z3 - z0;
786552f7358SJed Brown   const PetscScalar f_4       = x4 - x0;
787552f7358SJed Brown   const PetscScalar g_4       = y4 - y0;
788552f7358SJed Brown   const PetscScalar h_4       = z4 - z0;
789552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
790552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
791552f7358SJed Brown   const PetscScalar h_01      = z2 - z1 - z3 + z0;
792552f7358SJed Brown   const PetscScalar f_12      = x7 - x3 - x4 + x0;
793552f7358SJed Brown   const PetscScalar g_12      = y7 - y3 - y4 + y0;
794552f7358SJed Brown   const PetscScalar h_12      = z7 - z3 - z4 + z0;
795552f7358SJed Brown   const PetscScalar f_02      = x5 - x1 - x4 + x0;
796552f7358SJed Brown   const PetscScalar g_02      = y5 - y1 - y4 + y0;
797552f7358SJed Brown   const PetscScalar h_02      = z5 - z1 - z4 + z0;
798552f7358SJed Brown   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
799552f7358SJed Brown   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
800552f7358SJed Brown   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
80156044e6dSMatthew G. Knepley   const PetscScalar *ref;
80256044e6dSMatthew G. Knepley   PetscScalar       *real;
803552f7358SJed Brown   PetscErrorCode    ierr;
804552f7358SJed Brown 
805552f7358SJed Brown   PetscFunctionBegin;
80656044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
807552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
808552f7358SJed Brown   {
809552f7358SJed Brown     const PetscScalar p0 = ref[0];
810552f7358SJed Brown     const PetscScalar p1 = ref[1];
811552f7358SJed Brown     const PetscScalar p2 = ref[2];
812552f7358SJed Brown 
813552f7358SJed 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;
814552f7358SJed 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;
815552f7358SJed 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;
816552f7358SJed Brown   }
817552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
81856044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
819552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
820552f7358SJed Brown   PetscFunctionReturn(0);
821552f7358SJed Brown }
822552f7358SJed Brown 
823d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
824552f7358SJed Brown {
825552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
826552f7358SJed Brown   const PetscScalar x0        = vertices[0];
827552f7358SJed Brown   const PetscScalar y0        = vertices[1];
828552f7358SJed Brown   const PetscScalar z0        = vertices[2];
8297a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
8307a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
8317a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
832552f7358SJed Brown   const PetscScalar x2        = vertices[6];
833552f7358SJed Brown   const PetscScalar y2        = vertices[7];
834552f7358SJed Brown   const PetscScalar z2        = vertices[8];
8357a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
8367a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
8377a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
838552f7358SJed Brown   const PetscScalar x4        = vertices[12];
839552f7358SJed Brown   const PetscScalar y4        = vertices[13];
840552f7358SJed Brown   const PetscScalar z4        = vertices[14];
841552f7358SJed Brown   const PetscScalar x5        = vertices[15];
842552f7358SJed Brown   const PetscScalar y5        = vertices[16];
843552f7358SJed Brown   const PetscScalar z5        = vertices[17];
844552f7358SJed Brown   const PetscScalar x6        = vertices[18];
845552f7358SJed Brown   const PetscScalar y6        = vertices[19];
846552f7358SJed Brown   const PetscScalar z6        = vertices[20];
847552f7358SJed Brown   const PetscScalar x7        = vertices[21];
848552f7358SJed Brown   const PetscScalar y7        = vertices[22];
849552f7358SJed Brown   const PetscScalar z7        = vertices[23];
850552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
851552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
852552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
853552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
854552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
855552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
856552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
857552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
858552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
859552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
860552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
861552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
86256044e6dSMatthew G. Knepley   const PetscScalar *ref;
863552f7358SJed Brown   PetscErrorCode    ierr;
864552f7358SJed Brown 
865552f7358SJed Brown   PetscFunctionBegin;
86656044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
867552f7358SJed Brown   {
868552f7358SJed Brown     const PetscScalar x       = ref[0];
869552f7358SJed Brown     const PetscScalar y       = ref[1];
870552f7358SJed Brown     const PetscScalar z       = ref[2];
871552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
872da80777bSKarl Rupp     PetscScalar       values[9];
873da80777bSKarl Rupp 
874da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
875da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
876da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
877da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
878da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
879da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
880da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
881da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
882da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
8831aa26658SKarl Rupp 
88494ab13aaSBarry Smith     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
885552f7358SJed Brown   }
886552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
88756044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
88894ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88994ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
890552f7358SJed Brown   PetscFunctionReturn(0);
891552f7358SJed Brown }
892552f7358SJed Brown 
893a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
894a6dfd86eSKarl Rupp {
895fafc0619SMatthew G Knepley   DM             dmCoord;
896552f7358SJed Brown   SNES           snes;
897552f7358SJed Brown   KSP            ksp;
898552f7358SJed Brown   PC             pc;
899552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
900552f7358SJed Brown   Mat            J;
90156044e6dSMatthew G. Knepley   const PetscScalar *coords;
90256044e6dSMatthew G. Knepley   PetscScalar    *a;
903552f7358SJed Brown   PetscInt       p;
904552f7358SJed Brown   PetscErrorCode ierr;
905552f7358SJed Brown 
906552f7358SJed Brown   PetscFunctionBegin;
907552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
908fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
909552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
910552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
911552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
912552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
913c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
914552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
915552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
916552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
917552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
918552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
919552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
9200298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
9210298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
922552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
923552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
924552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
925552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
926552f7358SJed Brown 
92756044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
928552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
929552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
930a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
931552f7358SJed Brown     PetscScalar *xi;
932cb313848SJed Brown     PetscReal    xir[3];
933552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
934552f7358SJed Brown 
935552f7358SJed Brown     /* Can make this do all points at once */
9360298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
937ff1e0c32SBarry Smith     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
9380298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
939ff1e0c32SBarry Smith     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
9400298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
9410298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
942552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
943552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
944552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
945552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
946552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
947552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
948552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
949cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
950cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
951cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
952552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
953552f7358SJed Brown       a[p*ctx->dof+comp] =
954cb313848SJed Brown         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
9557a1931ceSMatthew G. Knepley         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
956cb313848SJed Brown         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
9577a1931ceSMatthew G. Knepley         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
958cb313848SJed Brown         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
959cb313848SJed Brown         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
960cb313848SJed Brown         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
961cb313848SJed Brown         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
962552f7358SJed Brown     }
963552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
9640298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
9650298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
966552f7358SJed Brown   }
967552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
96856044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
969552f7358SJed Brown 
970552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
971552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
972552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
973552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
974552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
975552f7358SJed Brown   PetscFunctionReturn(0);
976552f7358SJed Brown }
977552f7358SJed Brown 
9784267b1a3SMatthew G. Knepley /*@C
9794267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
9804267b1a3SMatthew G. Knepley 
981552f7358SJed Brown   Input Parameters:
982552f7358SJed Brown + ctx - The DMInterpolationInfo context
983552f7358SJed Brown . dm  - The DM
984552f7358SJed Brown - x   - The local vector containing the field to be interpolated
985552f7358SJed Brown 
986552f7358SJed Brown   Output Parameters:
987552f7358SJed Brown . v   - The vector containing the interpolated values
9884267b1a3SMatthew G. Knepley 
9894267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
9904267b1a3SMatthew G. Knepley 
9914267b1a3SMatthew G. Knepley   Level: beginner
9924267b1a3SMatthew G. Knepley 
9934267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
9944267b1a3SMatthew G. Knepley @*/
9950adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
9960adebc6cSBarry Smith {
997*680d18d5SMatthew G. Knepley   PetscInt       n;
998552f7358SJed Brown   PetscErrorCode ierr;
999552f7358SJed Brown 
1000552f7358SJed Brown   PetscFunctionBegin;
1001552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1002552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1003552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
1004552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1005ff1e0c32SBarry 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);
1006552f7358SJed Brown   if (n) {
1007*680d18d5SMatthew G. Knepley     PetscDS        ds;
1008*680d18d5SMatthew G. Knepley     DMPolytopeType ct;
1009*680d18d5SMatthew G. Knepley     PetscBool      done = PETSC_FALSE;
1010*680d18d5SMatthew G. Knepley 
1011*680d18d5SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1012*680d18d5SMatthew G. Knepley     if (ds) {
1013*680d18d5SMatthew G. Knepley       const PetscScalar *coords;
1014*680d18d5SMatthew G. Knepley       PetscScalar       *interpolant;
1015*680d18d5SMatthew G. Knepley       PetscInt           cdim, d, p, Nf, field, c = 0;
1016*680d18d5SMatthew G. Knepley 
1017*680d18d5SMatthew G. Knepley       ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1018*680d18d5SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
1019*680d18d5SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
1020*680d18d5SMatthew G. Knepley         PetscTabulation T;
1021*680d18d5SMatthew G. Knepley         PetscFE         fe;
1022*680d18d5SMatthew G. Knepley         PetscClassId    id;
1023*680d18d5SMatthew G. Knepley         PetscReal       xi[3];
1024*680d18d5SMatthew G. Knepley         PetscInt        Nc, f, fc;
1025*680d18d5SMatthew G. Knepley 
1026*680d18d5SMatthew G. Knepley         ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
1027*680d18d5SMatthew G. Knepley         ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
1028*680d18d5SMatthew G. Knepley         if (id != PETSCFE_CLASSID) break;
1029*680d18d5SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1030*680d18d5SMatthew G. Knepley         ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1031*680d18d5SMatthew G. Knepley         ierr = VecGetArrayWrite(v, &interpolant);CHKERRQ(ierr);
1032*680d18d5SMatthew G. Knepley         for (p = 0; p < ctx->n; ++p) {
1033*680d18d5SMatthew G. Knepley           PetscScalar *xa = NULL;
1034*680d18d5SMatthew G. Knepley           PetscReal    pcoords[3];
1035*680d18d5SMatthew G. Knepley 
1036*680d18d5SMatthew G. Knepley           for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
1037*680d18d5SMatthew G. Knepley           ierr = DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi);CHKERRQ(ierr);
1038*680d18d5SMatthew G. Knepley           ierr = DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1039*680d18d5SMatthew G. Knepley           ierr = PetscFECreateTabulation(fe, 1, 1, xi, 0, &T);CHKERRQ(ierr);
1040*680d18d5SMatthew G. Knepley           {
1041*680d18d5SMatthew G. Knepley             const PetscReal *basis = T->T[0];
1042*680d18d5SMatthew G. Knepley             const PetscInt   Nb    = T->Nb;
1043*680d18d5SMatthew G. Knepley             const PetscInt   Nc    = T->Nc;
1044*680d18d5SMatthew G. Knepley             for (fc = 0; fc < Nc; ++fc) {
1045*680d18d5SMatthew G. Knepley               interpolant[p*ctx->dof+c+fc] = 0.0;
1046*680d18d5SMatthew G. Knepley               for (f = 0; f < Nb; ++f) {
1047*680d18d5SMatthew G. Knepley                 interpolant[p*ctx->dof+c+fc] += xa[f]*basis[(0*Nb + f)*Nc + fc];
1048552f7358SJed Brown               }
1049*680d18d5SMatthew G. Knepley             }
1050*680d18d5SMatthew G. Knepley           }
1051*680d18d5SMatthew G. Knepley           ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
1052*680d18d5SMatthew G. Knepley           ierr = DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1053*680d18d5SMatthew G. Knepley         }
1054*680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayWrite(v, &interpolant);CHKERRQ(ierr);
1055*680d18d5SMatthew G. Knepley         ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1056*680d18d5SMatthew G. Knepley         c += Nc;
1057*680d18d5SMatthew G. Knepley       }
1058*680d18d5SMatthew G. Knepley       if (field == Nf) {
1059*680d18d5SMatthew G. Knepley         done = PETSC_TRUE;
1060*680d18d5SMatthew G. Knepley         if (c != ctx->dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %D != %D dof specified for interpolation", c, ctx->dof);
1061*680d18d5SMatthew G. Knepley       }
1062*680d18d5SMatthew G. Knepley     }
1063*680d18d5SMatthew G. Knepley     if (!done) {
1064*680d18d5SMatthew G. Knepley       /* TODO Check each cell individually */
1065*680d18d5SMatthew G. Knepley       ierr = DMPlexGetCellType(dm, ctx->cells[0], &ct);CHKERRQ(ierr);
1066*680d18d5SMatthew G. Knepley       switch (ct) {
1067*680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:      ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1068*680d18d5SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL: ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1069*680d18d5SMatthew G. Knepley         case DM_POLYTOPE_TETRAHEDRON:   ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1070*680d18d5SMatthew G. Knepley         case DM_POLYTOPE_HEXAHEDRON:    ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1071*680d18d5SMatthew G. Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support fpr cell type %s", DMPolytopeTypes[ct]);
1072*680d18d5SMatthew G. Knepley       }
1073*680d18d5SMatthew G. Knepley     }
1074552f7358SJed Brown   }
1075552f7358SJed Brown   PetscFunctionReturn(0);
1076552f7358SJed Brown }
1077552f7358SJed Brown 
10784267b1a3SMatthew G. Knepley /*@C
10794267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
10804267b1a3SMatthew G. Knepley 
10814267b1a3SMatthew G. Knepley   Collective on ctx
10824267b1a3SMatthew G. Knepley 
10834267b1a3SMatthew G. Knepley   Input Parameter:
10844267b1a3SMatthew G. Knepley . ctx - the context
10854267b1a3SMatthew G. Knepley 
10864267b1a3SMatthew G. Knepley   Level: beginner
10874267b1a3SMatthew G. Knepley 
10884267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
10894267b1a3SMatthew G. Knepley @*/
10900adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
10910adebc6cSBarry Smith {
1092552f7358SJed Brown   PetscErrorCode ierr;
1093552f7358SJed Brown 
1094552f7358SJed Brown   PetscFunctionBegin;
1095552f7358SJed Brown   PetscValidPointer(ctx, 2);
1096552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
1097552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
1098552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
1099552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
11000298fd71SBarry Smith   *ctx = NULL;
1101552f7358SJed Brown   PetscFunctionReturn(0);
1102552f7358SJed Brown }
1103cc0c4584SMatthew G. Knepley 
1104cc0c4584SMatthew G. Knepley /*@C
1105cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1106cc0c4584SMatthew G. Knepley 
1107cc0c4584SMatthew G. Knepley   Collective on SNES
1108cc0c4584SMatthew G. Knepley 
1109cc0c4584SMatthew G. Knepley   Input Parameters:
1110cc0c4584SMatthew G. Knepley + snes   - the SNES context
1111cc0c4584SMatthew G. Knepley . its    - iteration number
1112cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1113d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
1114cc0c4584SMatthew G. Knepley 
1115cc0c4584SMatthew G. Knepley   Notes:
1116cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1117cc0c4584SMatthew G. Knepley 
1118cc0c4584SMatthew G. Knepley   Level: intermediate
1119cc0c4584SMatthew G. Knepley 
1120cc0c4584SMatthew G. Knepley .seealso: SNESMonitorSet(), SNESMonitorDefault()
1121cc0c4584SMatthew G. Knepley @*/
1122d43b4f6eSBarry Smith PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1123cc0c4584SMatthew G. Knepley {
1124d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1125cc0c4584SMatthew G. Knepley   Vec                res;
1126cc0c4584SMatthew G. Knepley   DM                 dm;
1127cc0c4584SMatthew G. Knepley   PetscSection       s;
1128cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1129cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1130cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1131cc0c4584SMatthew G. Knepley   PetscErrorCode     ierr;
1132cc0c4584SMatthew G. Knepley 
1133cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11344d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
11359e5d0892SLisandro Dalcin   ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr);
1136cc0c4584SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
113792fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
1138cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
1139cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1140cc0c4584SMatthew G. Knepley   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
1141cc0c4584SMatthew G. Knepley   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
1142cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1143cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1144cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1145cc0c4584SMatthew G. Knepley 
1146cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1147cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1148cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1149cc0c4584SMatthew G. Knepley     }
1150cc0c4584SMatthew G. Knepley   }
1151cc0c4584SMatthew G. Knepley   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
1152b2566f29SBarry Smith   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1153d43b4f6eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1154cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1155cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
1156cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
1157cc0c4584SMatthew G. Knepley     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
1158cc0c4584SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
1159cc0c4584SMatthew G. Knepley   }
1160cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
1161cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1162d43b4f6eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1163cc0c4584SMatthew G. Knepley   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
1164cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1165cc0c4584SMatthew G. Knepley }
116624cdb843SMatthew G. Knepley 
116724cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
116824cdb843SMatthew G. Knepley 
116924cdb843SMatthew G. Knepley /*@
117024cdb843SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
117124cdb843SMatthew G. Knepley 
117224cdb843SMatthew G. Knepley   Input Parameters:
117324cdb843SMatthew G. Knepley + dm - The mesh
117424cdb843SMatthew G. Knepley . X  - Local solution
117524cdb843SMatthew G. Knepley - user - The user context
117624cdb843SMatthew G. Knepley 
117724cdb843SMatthew G. Knepley   Output Parameter:
117824cdb843SMatthew G. Knepley . F  - Local output vector
117924cdb843SMatthew G. Knepley 
118024cdb843SMatthew G. Knepley   Level: developer
118124cdb843SMatthew G. Knepley 
11827a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
118324cdb843SMatthew G. Knepley @*/
118424cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
118524cdb843SMatthew G. Knepley {
11866da023fcSToby Isaac   DM             plex;
1187083401c6SMatthew G. Knepley   IS             allcellIS;
1188083401c6SMatthew G. Knepley   PetscInt       Nds, s, depth;
118924cdb843SMatthew G. Knepley   PetscErrorCode ierr;
119024cdb843SMatthew G. Knepley 
119124cdb843SMatthew G. Knepley   PetscFunctionBegin;
1192083401c6SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
11936da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
11944a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1195083401c6SMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, &allcellIS);CHKERRQ(ierr);
1196083401c6SMatthew G. Knepley   if (!allcellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &allcellIS);CHKERRQ(ierr);}
1197083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1198083401c6SMatthew G. Knepley     PetscDS ds;
1199083401c6SMatthew G. Knepley     DMLabel label;
1200083401c6SMatthew G. Knepley     IS      cellIS;
1201083401c6SMatthew G. Knepley 
1202083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1203083401c6SMatthew G. Knepley     if (!label) {
1204083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1205083401c6SMatthew G. Knepley       cellIS = allcellIS;
1206083401c6SMatthew G. Knepley     } else {
1207083401c6SMatthew G. Knepley       IS pointIS;
1208083401c6SMatthew G. Knepley 
1209083401c6SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1210083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1211083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
12124a3e9fdbSToby Isaac     }
12134bee2e38SMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
12144a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1215083401c6SMatthew G. Knepley   }
1216083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
12179a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
121824cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
121924cdb843SMatthew G. Knepley }
122024cdb843SMatthew G. Knepley 
1221bdd6f66aSToby Isaac /*@
1222bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1223bdd6f66aSToby Isaac 
1224bdd6f66aSToby Isaac   Input Parameters:
1225bdd6f66aSToby Isaac + dm - The mesh
1226bdd6f66aSToby Isaac - user - The user context
1227bdd6f66aSToby Isaac 
1228bdd6f66aSToby Isaac   Output Parameter:
1229bdd6f66aSToby Isaac . X  - Local solution
1230bdd6f66aSToby Isaac 
1231bdd6f66aSToby Isaac   Level: developer
1232bdd6f66aSToby Isaac 
12337a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
1234bdd6f66aSToby Isaac @*/
1235bdd6f66aSToby Isaac PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1236bdd6f66aSToby Isaac {
1237bdd6f66aSToby Isaac   DM             plex;
1238bdd6f66aSToby Isaac   PetscErrorCode ierr;
1239bdd6f66aSToby Isaac 
1240bdd6f66aSToby Isaac   PetscFunctionBegin;
1241bdd6f66aSToby Isaac   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1242bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1243bdd6f66aSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1244bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1245bdd6f66aSToby Isaac }
1246bdd6f66aSToby Isaac 
12477a73cf09SMatthew G. Knepley /*@
12487a73cf09SMatthew G. Knepley   DMPlexComputeJacobianAction - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.
12497a73cf09SMatthew G. Knepley 
12507a73cf09SMatthew G. Knepley   Input Parameters:
12517a73cf09SMatthew G. Knepley + dm - The mesh
12527a73cf09SMatthew G. Knepley . cellIS -
12537a73cf09SMatthew G. Knepley . t  - The time
12547a73cf09SMatthew G. Knepley . X_tShift - The multiplier for the Jacobian with repsect to X_t
12557a73cf09SMatthew G. Knepley . X  - Local solution vector
12567a73cf09SMatthew G. Knepley . X_t  - Time-derivative of the local solution vector
12577a73cf09SMatthew G. Knepley . Y  - Local input vector
12587a73cf09SMatthew G. Knepley - user - The user context
12597a73cf09SMatthew G. Knepley 
12607a73cf09SMatthew G. Knepley   Output Parameter:
12617a73cf09SMatthew G. Knepley . Z - Local output vector
12627a73cf09SMatthew G. Knepley 
12637a73cf09SMatthew G. Knepley   Note:
12647a73cf09SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
12657a73cf09SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
12667a73cf09SMatthew G. Knepley 
12677a73cf09SMatthew G. Knepley   Level: developer
12687a73cf09SMatthew G. Knepley 
12697a73cf09SMatthew G. Knepley .seealso: FormFunctionLocal()
12707a73cf09SMatthew G. Knepley @*/
12717a73cf09SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
1272a925c78cSMatthew G. Knepley {
1273a925c78cSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1274a925c78cSMatthew G. Knepley   const char       *name  = "Jacobian";
12757a73cf09SMatthew G. Knepley   DM                dmAux, plex, plexAux = NULL;
1276a6e0b375SMatthew G. Knepley   DMEnclosureType   encAux;
1277c330f8ffSToby Isaac   Vec               A;
1278a925c78cSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
1279a925c78cSMatthew G. Knepley   PetscQuadrature   quad;
1280a925c78cSMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1281a925c78cSMatthew G. Knepley   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
12829044fa66SMatthew G. Knepley   PetscInt          Nf, fieldI, fieldJ;
12834d0b9603SSander Arens   PetscInt          totDim, totDimAux = 0;
12849044fa66SMatthew G. Knepley   const PetscInt   *cells;
12859044fa66SMatthew G. Knepley   PetscInt          cStart, cEnd, numCells, c;
128675b37a90SMatthew G. Knepley   PetscBool         hasDyn;
1287c330f8ffSToby Isaac   DMField           coordField;
1288a925c78cSMatthew G. Knepley   PetscErrorCode    ierr;
1289a925c78cSMatthew G. Knepley 
1290a925c78cSMatthew G. Knepley   PetscFunctionBegin;
1291a925c78cSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
12927a73cf09SMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
12937a73cf09SMatthew G. Knepley   if (!cellIS) {
12947a73cf09SMatthew G. Knepley     PetscInt depth;
12957a73cf09SMatthew G. Knepley 
12967a73cf09SMatthew G. Knepley     ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
12977a73cf09SMatthew G. Knepley     ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
12987a73cf09SMatthew G. Knepley     if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
12997a73cf09SMatthew G. Knepley   } else {
13007a73cf09SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) cellIS);CHKERRQ(ierr);
13017a73cf09SMatthew G. Knepley   }
130292fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
1303e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1304d28747ceSMatthew G. Knepley   ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr);
1305d28747ceSMatthew G. Knepley   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1306d28747ceSMatthew G. Knepley   ierr = DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);CHKERRQ(ierr);
1307d28747ceSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1308a925c78cSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1309a925c78cSMatthew G. Knepley   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
1310a925c78cSMatthew G. Knepley   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1311a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1312a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1313a925c78cSMatthew G. Knepley   if (dmAux) {
1314a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
13157a73cf09SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
131692fd8e1eSJed Brown     ierr = DMGetLocalSection(plexAux, &sectionAux);CHKERRQ(ierr);
1317a925c78cSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1318a925c78cSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1319a925c78cSMatthew G. Knepley   }
1320a925c78cSMatthew G. Knepley   ierr = VecSet(Z, 0.0);CHKERRQ(ierr);
1321a925c78cSMatthew G. Knepley   ierr = PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);CHKERRQ(ierr);
1322a925c78cSMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
13234a3e9fdbSToby Isaac   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1324a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
13259044fa66SMatthew G. Knepley     const PetscInt cell = cells ? cells[c] : c;
13269044fa66SMatthew G. Knepley     const PetscInt cind = c - cStart;
1327a925c78cSMatthew G. Knepley     PetscScalar   *x = NULL,  *x_t = NULL;
1328a925c78cSMatthew G. Knepley     PetscInt       i;
1329a925c78cSMatthew G. Knepley 
13309044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
13319044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
13329044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
1333a925c78cSMatthew G. Knepley     if (X_t) {
13349044fa66SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
13359044fa66SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
13369044fa66SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
1337a925c78cSMatthew G. Knepley     }
1338a925c78cSMatthew G. Knepley     if (dmAux) {
133944171101SMatthew G. Knepley       PetscInt subcell;
1340a6e0b375SMatthew G. Knepley       ierr = DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);CHKERRQ(ierr);
134144171101SMatthew G. Knepley       ierr = DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
13429044fa66SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
134344171101SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
1344a925c78cSMatthew G. Knepley     }
13459044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
13469044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
13479044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
1348a925c78cSMatthew G. Knepley   }
1349580bdb30SBarry Smith   ierr = PetscArrayzero(elemMat, numCells*totDim*totDim);CHKERRQ(ierr);
1350580bdb30SBarry Smith   if (hasDyn)  {ierr = PetscArrayzero(elemMatD, numCells*totDim*totDim);CHKERRQ(ierr);}
1351a925c78cSMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1352a925c78cSMatthew G. Knepley     PetscFE  fe;
1353c330f8ffSToby Isaac     PetscInt Nb;
1354a925c78cSMatthew G. Knepley     /* Conforming batches */
1355a925c78cSMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1356a925c78cSMatthew G. Knepley     /* Remainder */
1357c330f8ffSToby Isaac     PetscInt Nr, offset, Nq;
1358c330f8ffSToby Isaac     PetscQuadrature qGeom = NULL;
1359b7260050SToby Isaac     PetscInt    maxDegree;
1360c330f8ffSToby Isaac     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
1361a925c78cSMatthew G. Knepley 
1362a925c78cSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1363a925c78cSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1364a925c78cSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1365a925c78cSMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1366b7260050SToby Isaac     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
1367b7260050SToby Isaac     if (maxDegree <= 1) {ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);}
1368c330f8ffSToby Isaac     if (!qGeom) {
1369c330f8ffSToby Isaac       ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr);
1370c330f8ffSToby Isaac       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1371c330f8ffSToby Isaac     }
1372c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
13734a3e9fdbSToby Isaac     ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1374c330f8ffSToby Isaac     blockSize = Nb;
1375a925c78cSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1376a925c78cSMatthew G. Knepley     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1377a925c78cSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1378a925c78cSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1379a925c78cSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1380a925c78cSMatthew G. Knepley     offset    = numCells - Nr;
1381c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
1382c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1383a925c78cSMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
13844bee2e38SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
13854bee2e38SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
1386a925c78cSMatthew G. Knepley       if (hasDyn) {
13874bee2e38SMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
13884bee2e38SMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr);
1389a925c78cSMatthew G. Knepley       }
1390a925c78cSMatthew G. Knepley     }
1391c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1392c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
13934a3e9fdbSToby Isaac     ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1394c330f8ffSToby Isaac     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
1395a925c78cSMatthew G. Knepley   }
1396a925c78cSMatthew G. Knepley   if (hasDyn) {
13979044fa66SMatthew G. Knepley     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
1398a925c78cSMatthew G. Knepley   }
1399a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
14009044fa66SMatthew G. Knepley     const PetscInt     cell = cells ? cells[c] : c;
14019044fa66SMatthew G. Knepley     const PetscInt     cind = c - cStart;
1402a925c78cSMatthew G. Knepley     const PetscBLASInt M = totDim, one = 1;
1403a925c78cSMatthew G. Knepley     const PetscScalar  a = 1.0, b = 0.0;
1404a925c78cSMatthew G. Knepley 
14059044fa66SMatthew G. Knepley     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
1406a925c78cSMatthew G. Knepley     if (mesh->printFEM > 1) {
14079044fa66SMatthew G. Knepley       ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);
14089044fa66SMatthew G. Knepley       ierr = DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);CHKERRQ(ierr);
1409a925c78cSMatthew G. Knepley       ierr = DMPrintCellVector(c, "Z",  totDim, z);CHKERRQ(ierr);
1410a925c78cSMatthew G. Knepley     }
14119044fa66SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);CHKERRQ(ierr);
1412a925c78cSMatthew G. Knepley   }
1413a925c78cSMatthew G. Knepley   ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr);
1414a925c78cSMatthew G. Knepley   if (mesh->printFEM) {
1415ea13f565SStefano Zampini     ierr = PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n");CHKERRQ(ierr);
1416a8c5bd36SStefano Zampini     ierr = VecView(Z, NULL);CHKERRQ(ierr);
1417a925c78cSMatthew G. Knepley   }
14187a73cf09SMatthew G. Knepley   ierr = PetscFree(a);CHKERRQ(ierr);
14197a73cf09SMatthew G. Knepley   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
14207a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plexAux);CHKERRQ(ierr);
14217a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1422a925c78cSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1423a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1424a925c78cSMatthew G. Knepley }
1425a925c78cSMatthew G. Knepley 
142624cdb843SMatthew G. Knepley /*@
142724cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
142824cdb843SMatthew G. Knepley 
142924cdb843SMatthew G. Knepley   Input Parameters:
143024cdb843SMatthew G. Knepley + dm - The mesh
143124cdb843SMatthew G. Knepley . X  - Local input vector
143224cdb843SMatthew G. Knepley - user - The user context
143324cdb843SMatthew G. Knepley 
143424cdb843SMatthew G. Knepley   Output Parameter:
143524cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
143624cdb843SMatthew G. Knepley 
143724cdb843SMatthew G. Knepley   Note:
143824cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
143924cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
144024cdb843SMatthew G. Knepley 
144124cdb843SMatthew G. Knepley   Level: developer
144224cdb843SMatthew G. Knepley 
144324cdb843SMatthew G. Knepley .seealso: FormFunctionLocal()
144424cdb843SMatthew G. Knepley @*/
144524cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
144624cdb843SMatthew G. Knepley {
14476da023fcSToby Isaac   DM             plex;
1448083401c6SMatthew G. Knepley   IS             allcellIS;
1449f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
1450083401c6SMatthew G. Knepley   PetscInt       Nds, s, depth;
145124cdb843SMatthew G. Knepley   PetscErrorCode ierr;
145224cdb843SMatthew G. Knepley 
145324cdb843SMatthew G. Knepley   PetscFunctionBegin;
1454083401c6SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
14556da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
14564a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1457083401c6SMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, &allcellIS);CHKERRQ(ierr);
1458083401c6SMatthew G. Knepley   if (!allcellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &allcellIS);CHKERRQ(ierr);}
1459083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1460083401c6SMatthew G. Knepley     PetscDS ds;
1461083401c6SMatthew G. Knepley     DMLabel label;
1462083401c6SMatthew G. Knepley     IS      cellIS;
1463083401c6SMatthew G. Knepley 
1464083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1465083401c6SMatthew G. Knepley     if (!label) {
1466083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1467083401c6SMatthew G. Knepley       cellIS = allcellIS;
1468083401c6SMatthew G. Knepley     } else {
1469083401c6SMatthew G. Knepley       IS pointIS;
1470083401c6SMatthew G. Knepley 
1471083401c6SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1472083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1473083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1474083401c6SMatthew G. Knepley     }
1475083401c6SMatthew G. Knepley     if (!s) {
1476083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1477083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1478f04eb4edSMatthew G. Knepley       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
1479f04eb4edSMatthew G. Knepley       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1480083401c6SMatthew G. Knepley     }
14814a3e9fdbSToby Isaac     ierr = DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
14824a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1483083401c6SMatthew G. Knepley   }
1484083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
14859a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
148624cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
148724cdb843SMatthew G. Knepley }
14881878804aSMatthew G. Knepley 
148938cfc46eSPierre Jolivet /*
149038cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
149138cfc46eSPierre Jolivet 
149238cfc46eSPierre Jolivet    Input Parameters:
149338cfc46eSPierre Jolivet +     X - SNES linearization point
149438cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
149538cfc46eSPierre Jolivet 
149638cfc46eSPierre Jolivet    Output Parameter:
149738cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
149838cfc46eSPierre Jolivet 
149938cfc46eSPierre Jolivet    Level: intermediate
150038cfc46eSPierre Jolivet 
150138cfc46eSPierre Jolivet .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
150238cfc46eSPierre Jolivet */
150338cfc46eSPierre Jolivet static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
150438cfc46eSPierre Jolivet {
150538cfc46eSPierre Jolivet   SNES           snes;
150638cfc46eSPierre Jolivet   Mat            pJ;
150738cfc46eSPierre Jolivet   DM             ovldm,origdm;
150838cfc46eSPierre Jolivet   DMSNES         sdm;
150938cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM,Vec,void*);
151038cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
151138cfc46eSPierre Jolivet   void           *bctx,*jctx;
151238cfc46eSPierre Jolivet   PetscErrorCode ierr;
151338cfc46eSPierre Jolivet 
151438cfc46eSPierre Jolivet   PetscFunctionBegin;
151538cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr);
151638cfc46eSPierre Jolivet   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr);
151738cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr);
151838cfc46eSPierre Jolivet   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr);
151938cfc46eSPierre Jolivet   ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr);
152038cfc46eSPierre Jolivet   ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr);
152138cfc46eSPierre Jolivet   ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr);
152238cfc46eSPierre Jolivet   ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr);
152338cfc46eSPierre Jolivet   ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr);
152438cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr);
152538cfc46eSPierre Jolivet   if (!snes) {
152638cfc46eSPierre Jolivet     ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr);
152738cfc46eSPierre Jolivet     ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr);
152838cfc46eSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr);
152938cfc46eSPierre Jolivet     ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr);
153038cfc46eSPierre Jolivet   }
153138cfc46eSPierre Jolivet   ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr);
153238cfc46eSPierre Jolivet   ierr = VecLockReadPush(X);CHKERRQ(ierr);
153338cfc46eSPierre Jolivet   PetscStackPush("SNES user Jacobian function");
153438cfc46eSPierre Jolivet   ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr);
153538cfc46eSPierre Jolivet   PetscStackPop;
153638cfc46eSPierre Jolivet   ierr = VecLockReadPop(X);CHKERRQ(ierr);
153738cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
153838cfc46eSPierre Jolivet   {
153938cfc46eSPierre Jolivet     Mat locpJ;
154038cfc46eSPierre Jolivet 
154138cfc46eSPierre Jolivet     ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr);
154238cfc46eSPierre Jolivet     ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
154338cfc46eSPierre Jolivet   }
154438cfc46eSPierre Jolivet   PetscFunctionReturn(0);
154538cfc46eSPierre Jolivet }
154638cfc46eSPierre Jolivet 
1547a925c78cSMatthew G. Knepley /*@
15489f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
15499f520fc2SToby Isaac 
15509f520fc2SToby Isaac   Input Parameters:
15519f520fc2SToby Isaac + dm - The DM object
1552dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
15539f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
15549f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
15551a244344SSatish Balay 
15561a244344SSatish Balay   Level: developer
15579f520fc2SToby Isaac @*/
15589f520fc2SToby Isaac PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
15599f520fc2SToby Isaac {
15609f520fc2SToby Isaac   PetscErrorCode ierr;
15619f520fc2SToby Isaac 
15629f520fc2SToby Isaac   PetscFunctionBegin;
15639f520fc2SToby Isaac   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
15649f520fc2SToby Isaac   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
15659f520fc2SToby Isaac   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
156638cfc46eSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr);
15679f520fc2SToby Isaac   PetscFunctionReturn(0);
15689f520fc2SToby Isaac }
15699f520fc2SToby Isaac 
1570f017bcb6SMatthew G. Knepley /*@C
1571f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1572f017bcb6SMatthew G. Knepley 
1573f017bcb6SMatthew G. Knepley   Input Parameters:
1574f017bcb6SMatthew G. Knepley + snes - the SNES object
1575f017bcb6SMatthew G. Knepley . dm   - the DM
1576f2cacb80SMatthew G. Knepley . t    - the time
1577f017bcb6SMatthew G. Knepley . u    - a DM vector
1578f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1579f017bcb6SMatthew G. Knepley 
1580f3c94560SMatthew G. Knepley   Output Parameters:
1581f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1582f3c94560SMatthew G. Knepley 
15837f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
15847f96f943SMatthew G. Knepley 
1585f017bcb6SMatthew G. Knepley   Level: developer
1586f017bcb6SMatthew G. Knepley 
15877f96f943SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1588f017bcb6SMatthew G. Knepley @*/
1589f2cacb80SMatthew G. Knepley PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
15901878804aSMatthew G. Knepley {
1591f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1592f017bcb6SMatthew G. Knepley   void            **ectxs;
1593f3c94560SMatthew G. Knepley   PetscReal        *err;
15947f96f943SMatthew G. Knepley   MPI_Comm          comm;
15957f96f943SMatthew G. Knepley   PetscInt          Nf, f;
15961878804aSMatthew G. Knepley   PetscErrorCode    ierr;
15971878804aSMatthew G. Knepley 
15981878804aSMatthew G. Knepley   PetscFunctionBegin;
1599f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1600f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1601f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1602534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
16037f96f943SMatthew G. Knepley 
1604f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
16057f96f943SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
16067f96f943SMatthew G. Knepley 
1607f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1608f017bcb6SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1609083401c6SMatthew G. Knepley   ierr = PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr);
16107f96f943SMatthew G. Knepley   {
16117f96f943SMatthew G. Knepley     PetscInt Nds, s;
16127f96f943SMatthew G. Knepley 
1613083401c6SMatthew G. Knepley     ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1614083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
16157f96f943SMatthew G. Knepley       PetscDS         ds;
1616083401c6SMatthew G. Knepley       DMLabel         label;
1617083401c6SMatthew G. Knepley       IS              fieldIS;
16187f96f943SMatthew G. Knepley       const PetscInt *fields;
16197f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1620083401c6SMatthew G. Knepley 
1621083401c6SMatthew G. Knepley       ierr = DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
1622083401c6SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
1623083401c6SMatthew G. Knepley       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
1624083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1625083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
1626083401c6SMatthew G. Knepley         ierr = PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);CHKERRQ(ierr);
1627083401c6SMatthew G. Knepley       }
1628083401c6SMatthew G. Knepley       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
1629083401c6SMatthew G. Knepley     }
1630083401c6SMatthew G. Knepley   }
1631f017bcb6SMatthew G. Knepley   if (Nf > 1) {
1632f2cacb80SMatthew G. Knepley     ierr = DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);CHKERRQ(ierr);
1633f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1634f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1635f3c94560SMatthew 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);
16361878804aSMatthew G. Knepley       }
1637b878532fSMatthew G. Knepley     } else if (error) {
1638b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
16391878804aSMatthew G. Knepley     } else {
1640f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
1641f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1642f017bcb6SMatthew G. Knepley         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
1643f3c94560SMatthew G. Knepley         ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr);
16441878804aSMatthew G. Knepley       }
1645f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
1646f017bcb6SMatthew G. Knepley     }
1647f017bcb6SMatthew G. Knepley   } else {
1648f2cacb80SMatthew G. Knepley     ierr = DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);CHKERRQ(ierr);
1649f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1650f3c94560SMatthew G. Knepley       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1651b878532fSMatthew G. Knepley     } else if (error) {
1652b878532fSMatthew G. Knepley       error[0] = err[0];
1653f017bcb6SMatthew G. Knepley     } else {
1654f3c94560SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);CHKERRQ(ierr);
1655f017bcb6SMatthew G. Knepley     }
1656f017bcb6SMatthew G. Knepley   }
1657f3c94560SMatthew G. Knepley   ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr);
1658f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1659f017bcb6SMatthew G. Knepley }
1660f017bcb6SMatthew G. Knepley 
1661f017bcb6SMatthew G. Knepley /*@C
1662f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1663f017bcb6SMatthew G. Knepley 
1664f017bcb6SMatthew G. Knepley   Input Parameters:
1665f017bcb6SMatthew G. Knepley + snes - the SNES object
1666f017bcb6SMatthew G. Knepley . dm   - the DM
1667f017bcb6SMatthew G. Knepley . u    - a DM vector
1668f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1669f017bcb6SMatthew G. Knepley 
1670f3c94560SMatthew G. Knepley   Output Parameters:
1671f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1672f3c94560SMatthew G. Knepley 
1673f017bcb6SMatthew G. Knepley   Level: developer
1674f017bcb6SMatthew G. Knepley 
1675f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1676f017bcb6SMatthew G. Knepley @*/
1677f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1678f017bcb6SMatthew G. Knepley {
1679f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1680f017bcb6SMatthew G. Knepley   Vec            r;
1681f017bcb6SMatthew G. Knepley   PetscReal      res;
1682f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1683f017bcb6SMatthew G. Knepley 
1684f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1685f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1686f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1687f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1688534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
1689f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1690f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1691f017bcb6SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
16921878804aSMatthew G. Knepley   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
16931878804aSMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
1694f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
1695f017bcb6SMatthew G. Knepley     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1696b878532fSMatthew G. Knepley   } else if (residual) {
1697b878532fSMatthew G. Knepley     *residual = res;
1698f017bcb6SMatthew G. Knepley   } else {
1699f017bcb6SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
17001878804aSMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
17011878804aSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
1702685405a1SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
1703685405a1SBarry Smith     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
1704f017bcb6SMatthew G. Knepley   }
1705f017bcb6SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
1706f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1707f017bcb6SMatthew G. Knepley }
1708f017bcb6SMatthew G. Knepley 
1709f017bcb6SMatthew G. Knepley /*@C
1710f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1711f017bcb6SMatthew G. Knepley 
1712f017bcb6SMatthew G. Knepley   Input Parameters:
1713f017bcb6SMatthew G. Knepley + snes - the SNES object
1714f017bcb6SMatthew G. Knepley . dm   - the DM
1715f017bcb6SMatthew G. Knepley . u    - a DM vector
1716f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1717f017bcb6SMatthew G. Knepley 
1718f3c94560SMatthew G. Knepley   Output Parameters:
1719f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1720f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1721f3c94560SMatthew G. Knepley 
1722f017bcb6SMatthew G. Knepley   Level: developer
1723f017bcb6SMatthew G. Knepley 
1724f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1725f017bcb6SMatthew G. Knepley @*/
1726f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1727f017bcb6SMatthew G. Knepley {
1728f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1729f017bcb6SMatthew G. Knepley   PetscDS        ds;
1730f017bcb6SMatthew G. Knepley   Mat            J, M;
1731f017bcb6SMatthew G. Knepley   MatNullSpace   nullspace;
1732f3c94560SMatthew G. Knepley   PetscReal      slope, intercept;
17337f96f943SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1734f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1735f017bcb6SMatthew G. Knepley 
1736f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1737f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1738f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1739f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1740534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1741534a8f05SLisandro Dalcin   if (convRate) PetscValidRealPointer(convRate, 5);
1742f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1743f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1744f017bcb6SMatthew G. Knepley   /* Create and view matrices */
1745f017bcb6SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
17467f96f943SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1747f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1748f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1749282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
1750282e7bb4SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
1751282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
1752f017bcb6SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
1753282e7bb4SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
1754282e7bb4SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
1755282e7bb4SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
1756282e7bb4SMatthew G. Knepley   } else {
1757282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
1758282e7bb4SMatthew G. Knepley   }
1759f017bcb6SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
1760282e7bb4SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
1761282e7bb4SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
1762f017bcb6SMatthew G. Knepley   /* Check nullspace */
1763f017bcb6SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
1764f017bcb6SMatthew G. Knepley   if (nullspace) {
17651878804aSMatthew G. Knepley     PetscBool isNull;
1766f017bcb6SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
1767f017bcb6SMatthew G. Knepley     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
17681878804aSMatthew G. Knepley   }
1769f017bcb6SMatthew G. Knepley   /* Taylor test */
1770f017bcb6SMatthew G. Knepley   {
1771f017bcb6SMatthew G. Knepley     PetscRandom rand;
1772f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1773f3c94560SMatthew G. Knepley     PetscReal   h;
1774f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1775f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1776f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1777f017bcb6SMatthew G. Knepley 
1778f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
1779f017bcb6SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
1780f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
1781f017bcb6SMatthew G. Knepley     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
1782f017bcb6SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
1783f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
1784f017bcb6SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
1785f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
1786f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1787f017bcb6SMatthew G. Knepley     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1788f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1789f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1790f017bcb6SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
1791f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
1792f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
1793f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1794f017bcb6SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
1795f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1796f017bcb6SMatthew G. Knepley       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
1797f017bcb6SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
1798f017bcb6SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
1799f017bcb6SMatthew G. Knepley 
1800f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1801f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1802f017bcb6SMatthew G. Knepley     }
1803f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
1804f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
1805f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
18061878804aSMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
1807f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
1808f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1809f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1810f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
1811f017bcb6SMatthew G. Knepley     }
1812f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
1813f017bcb6SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
1814f017bcb6SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
1815f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1816f017bcb6SMatthew G. Knepley     if (tol >= 0) {
1817b878532fSMatthew 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);
1818b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1819b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1820b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1821f017bcb6SMatthew G. Knepley     } else {
1822f3c94560SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
1823f017bcb6SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
1824f017bcb6SMatthew G. Knepley     }
1825f017bcb6SMatthew G. Knepley   }
18261878804aSMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
18271878804aSMatthew G. Knepley   PetscFunctionReturn(0);
18281878804aSMatthew G. Knepley }
18291878804aSMatthew G. Knepley 
18307f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1831f017bcb6SMatthew G. Knepley {
1832f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1833f017bcb6SMatthew G. Knepley 
1834f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1835f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);CHKERRQ(ierr);
1836f3c94560SMatthew G. Knepley   ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr);
1837f3c94560SMatthew G. Knepley   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr);
1838f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1839f017bcb6SMatthew G. Knepley }
1840f017bcb6SMatthew G. Knepley 
1841bee9a294SMatthew G. Knepley /*@C
1842bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1843bee9a294SMatthew G. Knepley 
1844bee9a294SMatthew G. Knepley   Input Parameters:
1845bee9a294SMatthew G. Knepley + snes - the SNES object
18467f96f943SMatthew G. Knepley - u    - representative SNES vector
18477f96f943SMatthew G. Knepley 
18487f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1849bee9a294SMatthew G. Knepley 
1850bee9a294SMatthew G. Knepley   Level: developer
1851bee9a294SMatthew G. Knepley @*/
18527f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
18531878804aSMatthew G. Knepley {
18541878804aSMatthew G. Knepley   DM             dm;
18551878804aSMatthew G. Knepley   Vec            sol;
18561878804aSMatthew G. Knepley   PetscBool      check;
18571878804aSMatthew G. Knepley   PetscErrorCode ierr;
18581878804aSMatthew G. Knepley 
18591878804aSMatthew G. Knepley   PetscFunctionBegin;
1860c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
18611878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
18621878804aSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
18631878804aSMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
18641878804aSMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
18657f96f943SMatthew G. Knepley   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
18661878804aSMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1867552f7358SJed Brown   PetscFunctionReturn(0);
1868552f7358SJed Brown }
1869