xref: /petsc/src/snes/utils/dmplexsnes.c (revision a1e447459e008a1aada86f00c8024fda75255070)
1552f7358SJed Brown #include <petscdmplex.h> /*I "petscdmplex.h" I*/
2552f7358SJed Brown #include <petscsnes.h>      /*I "petscsnes.h" I*/
3afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>
4552f7358SJed Brown 
5552f7358SJed Brown #undef __FUNCT__
6552f7358SJed Brown #define __FUNCT__ "DMInterpolationCreate"
70adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
80adebc6cSBarry Smith {
9552f7358SJed Brown   PetscErrorCode ierr;
10552f7358SJed Brown 
11552f7358SJed Brown   PetscFunctionBegin;
12552f7358SJed Brown   PetscValidPointer(ctx, 2);
13552f7358SJed Brown   ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr);
141aa26658SKarl Rupp 
15552f7358SJed Brown   (*ctx)->comm   = comm;
16552f7358SJed Brown   (*ctx)->dim    = -1;
17552f7358SJed Brown   (*ctx)->nInput = 0;
180298fd71SBarry Smith   (*ctx)->points = NULL;
190298fd71SBarry Smith   (*ctx)->cells  = NULL;
20552f7358SJed Brown   (*ctx)->n      = -1;
210298fd71SBarry Smith   (*ctx)->coords = NULL;
22552f7358SJed Brown   PetscFunctionReturn(0);
23552f7358SJed Brown }
24552f7358SJed Brown 
25552f7358SJed Brown #undef __FUNCT__
26552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDim"
270adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
280adebc6cSBarry Smith {
29552f7358SJed Brown   PetscFunctionBegin;
300adebc6cSBarry Smith   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);
31552f7358SJed Brown   ctx->dim = dim;
32552f7358SJed Brown   PetscFunctionReturn(0);
33552f7358SJed Brown }
34552f7358SJed Brown 
35552f7358SJed Brown #undef __FUNCT__
36552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDim"
370adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
380adebc6cSBarry Smith {
39552f7358SJed Brown   PetscFunctionBegin;
40552f7358SJed Brown   PetscValidIntPointer(dim, 2);
41552f7358SJed Brown   *dim = ctx->dim;
42552f7358SJed Brown   PetscFunctionReturn(0);
43552f7358SJed Brown }
44552f7358SJed Brown 
45552f7358SJed Brown #undef __FUNCT__
46552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDof"
470adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
480adebc6cSBarry Smith {
49552f7358SJed Brown   PetscFunctionBegin;
500adebc6cSBarry Smith   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);
51552f7358SJed Brown   ctx->dof = dof;
52552f7358SJed Brown   PetscFunctionReturn(0);
53552f7358SJed Brown }
54552f7358SJed Brown 
55552f7358SJed Brown #undef __FUNCT__
56552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDof"
570adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
580adebc6cSBarry Smith {
59552f7358SJed Brown   PetscFunctionBegin;
60552f7358SJed Brown   PetscValidIntPointer(dof, 2);
61552f7358SJed Brown   *dof = ctx->dof;
62552f7358SJed Brown   PetscFunctionReturn(0);
63552f7358SJed Brown }
64552f7358SJed Brown 
65552f7358SJed Brown #undef __FUNCT__
66552f7358SJed Brown #define __FUNCT__ "DMInterpolationAddPoints"
670adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
680adebc6cSBarry Smith {
69552f7358SJed Brown   PetscErrorCode ierr;
70552f7358SJed Brown 
71552f7358SJed Brown   PetscFunctionBegin;
720adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
730adebc6cSBarry Smith   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
74552f7358SJed Brown   ctx->nInput = n;
751aa26658SKarl Rupp 
76552f7358SJed Brown   ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr);
77552f7358SJed Brown   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
78552f7358SJed Brown   PetscFunctionReturn(0);
79552f7358SJed Brown }
80552f7358SJed Brown 
81552f7358SJed Brown #undef __FUNCT__
82552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetUp"
830adebc6cSBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
840adebc6cSBarry Smith {
85552f7358SJed Brown   MPI_Comm       comm = ctx->comm;
86552f7358SJed Brown   PetscScalar    *a;
87552f7358SJed Brown   PetscInt       p, q, i;
88552f7358SJed Brown   PetscMPIInt    rank, size;
89552f7358SJed Brown   PetscErrorCode ierr;
90552f7358SJed Brown   Vec            pointVec;
91552f7358SJed Brown   IS             cellIS;
92552f7358SJed Brown   PetscLayout    layout;
93552f7358SJed Brown   PetscReal      *globalPoints;
94cb313848SJed Brown   PetscScalar    *globalPointsScalar;
95552f7358SJed Brown   const PetscInt *ranges;
96552f7358SJed Brown   PetscMPIInt    *counts, *displs;
97552f7358SJed Brown   const PetscInt *foundCells;
98552f7358SJed Brown   PetscMPIInt    *foundProcs, *globalProcs;
9919436ca2SJed Brown   PetscInt       n, N;
100552f7358SJed Brown 
10119436ca2SJed Brown   PetscFunctionBegin;
10219436ca2SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10319436ca2SJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10419436ca2SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1050adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
10619436ca2SJed Brown   /* Locate points */
10719436ca2SJed Brown   n = ctx->nInput;
108552f7358SJed Brown   if (!redundantPoints) {
109552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
110552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
111552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
112552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
113552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
114552f7358SJed Brown     /* Communicate all points to all processes */
115552f7358SJed Brown     ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr);
116552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
117552f7358SJed Brown     for (p = 0; p < size; ++p) {
118552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
119552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
120552f7358SJed Brown     }
121552f7358SJed Brown     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
122552f7358SJed Brown   } else {
123552f7358SJed Brown     N = n;
124f5af7f23SKarl Rupp 
125552f7358SJed Brown     globalPoints = ctx->points;
126552f7358SJed Brown   }
127552f7358SJed Brown #if 0
128552f7358SJed Brown   ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
12919436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
130552f7358SJed Brown #else
131cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
132cb313848SJed Brown   ierr = PetscMalloc(N*sizeof(PetscScalar),&globalPointsScalar);CHKERRQ(ierr);
133cb313848SJed Brown   for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i];
134cb313848SJed Brown #else
135cb313848SJed Brown   globalPointsScalar = globalPoints;
136cb313848SJed Brown #endif
13704706141SMatthew G Knepley   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
138552f7358SJed Brown   ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
139552f7358SJed Brown   ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr);
140552f7358SJed Brown   ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr);
141552f7358SJed Brown #endif
142552f7358SJed Brown   for (p = 0; p < N; ++p) {
1431aa26658SKarl Rupp     if (foundCells[p] >= 0) foundProcs[p] = rank;
1441aa26658SKarl Rupp     else foundProcs[p] = size;
145552f7358SJed Brown   }
146552f7358SJed Brown   /* Let the lowest rank process own each point */
147efab3cc2SJed Brown   ierr   = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
148552f7358SJed Brown   ctx->n = 0;
149552f7358SJed Brown   for (p = 0; p < N; ++p) {
1500adebc6cSBarry Smith     if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, globalPoints[p*ctx->dim+0], ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0, ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0);
1511aa26658SKarl Rupp     else if (globalProcs[p] == rank) ctx->n++;
152552f7358SJed Brown   }
153552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
154552f7358SJed Brown   ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr);
155552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
156552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
157552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
158552f7358SJed Brown   ierr = VecSetFromOptions(ctx->coords);CHKERRQ(ierr);
159552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
160552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
161552f7358SJed Brown     if (globalProcs[p] == rank) {
162552f7358SJed Brown       PetscInt d;
163552f7358SJed Brown 
1641aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
165552f7358SJed Brown       ctx->cells[q++] = foundCells[p];
166552f7358SJed Brown     }
167552f7358SJed Brown   }
168552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
169552f7358SJed Brown #if 0
170552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
171552f7358SJed Brown #else
172552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
173552f7358SJed Brown   ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr);
174552f7358SJed Brown   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
175552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
176552f7358SJed Brown #endif
177cb313848SJed Brown   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
178552f7358SJed Brown   if (!redundantPoints) {
179552f7358SJed Brown     ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);
180552f7358SJed Brown     ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
181552f7358SJed Brown   }
182552f7358SJed Brown   PetscFunctionReturn(0);
183552f7358SJed Brown }
184552f7358SJed Brown 
185552f7358SJed Brown #undef __FUNCT__
186552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetCoordinates"
1870adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
1880adebc6cSBarry Smith {
189552f7358SJed Brown   PetscFunctionBegin;
190552f7358SJed Brown   PetscValidPointer(coordinates, 2);
1910adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
192552f7358SJed Brown   *coordinates = ctx->coords;
193552f7358SJed Brown   PetscFunctionReturn(0);
194552f7358SJed Brown }
195552f7358SJed Brown 
196552f7358SJed Brown #undef __FUNCT__
197552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetVector"
1980adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
1990adebc6cSBarry Smith {
200552f7358SJed Brown   PetscErrorCode ierr;
201552f7358SJed Brown 
202552f7358SJed Brown   PetscFunctionBegin;
203552f7358SJed Brown   PetscValidPointer(v, 2);
2040adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
205552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
206552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
207552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
208552f7358SJed Brown   ierr = VecSetFromOptions(*v);CHKERRQ(ierr);
209552f7358SJed Brown   PetscFunctionReturn(0);
210552f7358SJed Brown }
211552f7358SJed Brown 
212552f7358SJed Brown #undef __FUNCT__
213552f7358SJed Brown #define __FUNCT__ "DMInterpolationRestoreVector"
2140adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
2150adebc6cSBarry Smith {
216552f7358SJed Brown   PetscErrorCode ierr;
217552f7358SJed Brown 
218552f7358SJed Brown   PetscFunctionBegin;
219552f7358SJed Brown   PetscValidPointer(v, 2);
2200adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
221552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
222552f7358SJed Brown   PetscFunctionReturn(0);
223552f7358SJed Brown }
224552f7358SJed Brown 
225552f7358SJed Brown #undef __FUNCT__
226552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Simplex_Private"
227a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
228a6dfd86eSKarl Rupp {
229552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
230552f7358SJed Brown   PetscScalar    *a, *coords;
231552f7358SJed Brown   PetscInt       p;
232552f7358SJed Brown   PetscErrorCode ierr;
233552f7358SJed Brown 
234552f7358SJed Brown   PetscFunctionBegin;
235552f7358SJed Brown   ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr);
236552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
237552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
238552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
239552f7358SJed Brown     PetscInt     c = ctx->cells[p];
240*a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
241552f7358SJed Brown     PetscReal    xi[4];
242552f7358SJed Brown     PetscInt     d, f, comp;
243552f7358SJed Brown 
244552f7358SJed Brown     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
245552f7358SJed Brown     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
2460298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
2471aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
2481aa26658SKarl Rupp 
249552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
250552f7358SJed Brown       xi[d] = 0.0;
2511aa26658SKarl 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]);
2521aa26658SKarl 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];
253552f7358SJed Brown     }
2540298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
255552f7358SJed Brown   }
256552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
257552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
258552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
259552f7358SJed Brown   PetscFunctionReturn(0);
260552f7358SJed Brown }
261552f7358SJed Brown 
262552f7358SJed Brown #undef __FUNCT__
263552f7358SJed Brown #define __FUNCT__ "QuadMap_Private"
2645820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
265552f7358SJed Brown {
266552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
267552f7358SJed Brown   const PetscScalar x0        = vertices[0];
268552f7358SJed Brown   const PetscScalar y0        = vertices[1];
269552f7358SJed Brown   const PetscScalar x1        = vertices[2];
270552f7358SJed Brown   const PetscScalar y1        = vertices[3];
271552f7358SJed Brown   const PetscScalar x2        = vertices[4];
272552f7358SJed Brown   const PetscScalar y2        = vertices[5];
273552f7358SJed Brown   const PetscScalar x3        = vertices[6];
274552f7358SJed Brown   const PetscScalar y3        = vertices[7];
275552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
276552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
277552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
278552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
279552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
280552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
281552f7358SJed Brown   PetscScalar       *ref, *real;
282552f7358SJed Brown   PetscErrorCode    ierr;
283552f7358SJed Brown 
284552f7358SJed Brown   PetscFunctionBegin;
285552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
286552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
287552f7358SJed Brown   {
288552f7358SJed Brown     const PetscScalar p0 = ref[0];
289552f7358SJed Brown     const PetscScalar p1 = ref[1];
290552f7358SJed Brown 
291552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
292552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
293552f7358SJed Brown   }
294552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
295552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
296552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
297552f7358SJed Brown   PetscFunctionReturn(0);
298552f7358SJed Brown }
299552f7358SJed Brown 
300552f7358SJed Brown #undef __FUNCT__
301552f7358SJed Brown #define __FUNCT__ "QuadJacobian_Private"
3025820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
303552f7358SJed Brown {
304552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
305552f7358SJed Brown   const PetscScalar x0        = vertices[0];
306552f7358SJed Brown   const PetscScalar y0        = vertices[1];
307552f7358SJed Brown   const PetscScalar x1        = vertices[2];
308552f7358SJed Brown   const PetscScalar y1        = vertices[3];
309552f7358SJed Brown   const PetscScalar x2        = vertices[4];
310552f7358SJed Brown   const PetscScalar y2        = vertices[5];
311552f7358SJed Brown   const PetscScalar x3        = vertices[6];
312552f7358SJed Brown   const PetscScalar y3        = vertices[7];
313552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
314552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
315552f7358SJed Brown   PetscScalar       *ref;
316552f7358SJed Brown   PetscErrorCode    ierr;
317552f7358SJed Brown 
318552f7358SJed Brown   PetscFunctionBegin;
319552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
320552f7358SJed Brown   {
321552f7358SJed Brown     const PetscScalar x       = ref[0];
322552f7358SJed Brown     const PetscScalar y       = ref[1];
323552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
324da80777bSKarl Rupp     PetscScalar       values[4];
325da80777bSKarl Rupp 
326da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
327da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
328552f7358SJed Brown     ierr      = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
329552f7358SJed Brown   }
330552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
331552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
332552f7358SJed Brown   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333552f7358SJed Brown   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334552f7358SJed Brown   PetscFunctionReturn(0);
335552f7358SJed Brown }
336552f7358SJed Brown 
337552f7358SJed Brown #undef __FUNCT__
338552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Quad_Private"
339a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
340a6dfd86eSKarl Rupp {
341fafc0619SMatthew G Knepley   DM             dmCoord;
342552f7358SJed Brown   SNES           snes;
343552f7358SJed Brown   KSP            ksp;
344552f7358SJed Brown   PC             pc;
345552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
346552f7358SJed Brown   Mat            J;
347552f7358SJed Brown   PetscScalar    *a, *coords;
348552f7358SJed Brown   PetscInt       p;
349552f7358SJed Brown   PetscErrorCode ierr;
350552f7358SJed Brown 
351552f7358SJed Brown   PetscFunctionBegin;
352552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
353fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
354552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
355552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
356552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
357552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
358552f7358SJed Brown   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
359552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
360552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
361552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
362552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
363552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
364552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
3650298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
3660298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
367552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
368552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
369552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
370552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
371552f7358SJed Brown 
372552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
373552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
374552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
375*a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
376552f7358SJed Brown     PetscScalar *xi;
377cb313848SJed Brown     PetscReal    xir[2];
378552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
379552f7358SJed Brown 
380552f7358SJed Brown     /* Can make this do all points at once */
3810298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
3820adebc6cSBarry Smith     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
3830298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
3840adebc6cSBarry Smith     if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);
3850298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
3860298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
387552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
388552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
389552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
390552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
391552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
392552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
393cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
394cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
3951aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*ctx->dof+comp]*xir[0]*(1 - xir[1]) + x[2*ctx->dof+comp]*xir[0]*xir[1] + x[3*ctx->dof+comp]*(1 - xir[0])*xir[1];
3961aa26658SKarl Rupp 
397552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
3980298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
3990298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
400552f7358SJed Brown   }
401552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
402552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
403552f7358SJed Brown 
404552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
405552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
406552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
407552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
408552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
409552f7358SJed Brown   PetscFunctionReturn(0);
410552f7358SJed Brown }
411552f7358SJed Brown 
412552f7358SJed Brown #undef __FUNCT__
413552f7358SJed Brown #define __FUNCT__ "HexMap_Private"
4145820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
415552f7358SJed Brown {
416552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
417552f7358SJed Brown   const PetscScalar x0        = vertices[0];
418552f7358SJed Brown   const PetscScalar y0        = vertices[1];
419552f7358SJed Brown   const PetscScalar z0        = vertices[2];
420552f7358SJed Brown   const PetscScalar x1        = vertices[3];
421552f7358SJed Brown   const PetscScalar y1        = vertices[4];
422552f7358SJed Brown   const PetscScalar z1        = vertices[5];
423552f7358SJed Brown   const PetscScalar x2        = vertices[6];
424552f7358SJed Brown   const PetscScalar y2        = vertices[7];
425552f7358SJed Brown   const PetscScalar z2        = vertices[8];
426552f7358SJed Brown   const PetscScalar x3        = vertices[9];
427552f7358SJed Brown   const PetscScalar y3        = vertices[10];
428552f7358SJed Brown   const PetscScalar z3        = vertices[11];
429552f7358SJed Brown   const PetscScalar x4        = vertices[12];
430552f7358SJed Brown   const PetscScalar y4        = vertices[13];
431552f7358SJed Brown   const PetscScalar z4        = vertices[14];
432552f7358SJed Brown   const PetscScalar x5        = vertices[15];
433552f7358SJed Brown   const PetscScalar y5        = vertices[16];
434552f7358SJed Brown   const PetscScalar z5        = vertices[17];
435552f7358SJed Brown   const PetscScalar x6        = vertices[18];
436552f7358SJed Brown   const PetscScalar y6        = vertices[19];
437552f7358SJed Brown   const PetscScalar z6        = vertices[20];
438552f7358SJed Brown   const PetscScalar x7        = vertices[21];
439552f7358SJed Brown   const PetscScalar y7        = vertices[22];
440552f7358SJed Brown   const PetscScalar z7        = vertices[23];
441552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
442552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
443552f7358SJed Brown   const PetscScalar h_1       = z1 - z0;
444552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
445552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
446552f7358SJed Brown   const PetscScalar h_3       = z3 - z0;
447552f7358SJed Brown   const PetscScalar f_4       = x4 - x0;
448552f7358SJed Brown   const PetscScalar g_4       = y4 - y0;
449552f7358SJed Brown   const PetscScalar h_4       = z4 - z0;
450552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
451552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
452552f7358SJed Brown   const PetscScalar h_01      = z2 - z1 - z3 + z0;
453552f7358SJed Brown   const PetscScalar f_12      = x7 - x3 - x4 + x0;
454552f7358SJed Brown   const PetscScalar g_12      = y7 - y3 - y4 + y0;
455552f7358SJed Brown   const PetscScalar h_12      = z7 - z3 - z4 + z0;
456552f7358SJed Brown   const PetscScalar f_02      = x5 - x1 - x4 + x0;
457552f7358SJed Brown   const PetscScalar g_02      = y5 - y1 - y4 + y0;
458552f7358SJed Brown   const PetscScalar h_02      = z5 - z1 - z4 + z0;
459552f7358SJed Brown   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
460552f7358SJed Brown   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
461552f7358SJed Brown   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
462552f7358SJed Brown   PetscScalar       *ref, *real;
463552f7358SJed Brown   PetscErrorCode    ierr;
464552f7358SJed Brown 
465552f7358SJed Brown   PetscFunctionBegin;
466552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
467552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
468552f7358SJed Brown   {
469552f7358SJed Brown     const PetscScalar p0 = ref[0];
470552f7358SJed Brown     const PetscScalar p1 = ref[1];
471552f7358SJed Brown     const PetscScalar p2 = ref[2];
472552f7358SJed Brown 
473552f7358SJed 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;
474552f7358SJed 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;
475552f7358SJed 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;
476552f7358SJed Brown   }
477552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
478552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
479552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
480552f7358SJed Brown   PetscFunctionReturn(0);
481552f7358SJed Brown }
482552f7358SJed Brown 
483552f7358SJed Brown #undef __FUNCT__
484552f7358SJed Brown #define __FUNCT__ "HexJacobian_Private"
4855820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
486552f7358SJed Brown {
487552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
488552f7358SJed Brown   const PetscScalar x0        = vertices[0];
489552f7358SJed Brown   const PetscScalar y0        = vertices[1];
490552f7358SJed Brown   const PetscScalar z0        = vertices[2];
491552f7358SJed Brown   const PetscScalar x1        = vertices[3];
492552f7358SJed Brown   const PetscScalar y1        = vertices[4];
493552f7358SJed Brown   const PetscScalar z1        = vertices[5];
494552f7358SJed Brown   const PetscScalar x2        = vertices[6];
495552f7358SJed Brown   const PetscScalar y2        = vertices[7];
496552f7358SJed Brown   const PetscScalar z2        = vertices[8];
497552f7358SJed Brown   const PetscScalar x3        = vertices[9];
498552f7358SJed Brown   const PetscScalar y3        = vertices[10];
499552f7358SJed Brown   const PetscScalar z3        = vertices[11];
500552f7358SJed Brown   const PetscScalar x4        = vertices[12];
501552f7358SJed Brown   const PetscScalar y4        = vertices[13];
502552f7358SJed Brown   const PetscScalar z4        = vertices[14];
503552f7358SJed Brown   const PetscScalar x5        = vertices[15];
504552f7358SJed Brown   const PetscScalar y5        = vertices[16];
505552f7358SJed Brown   const PetscScalar z5        = vertices[17];
506552f7358SJed Brown   const PetscScalar x6        = vertices[18];
507552f7358SJed Brown   const PetscScalar y6        = vertices[19];
508552f7358SJed Brown   const PetscScalar z6        = vertices[20];
509552f7358SJed Brown   const PetscScalar x7        = vertices[21];
510552f7358SJed Brown   const PetscScalar y7        = vertices[22];
511552f7358SJed Brown   const PetscScalar z7        = vertices[23];
512552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
513552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
514552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
515552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
516552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
517552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
518552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
519552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
520552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
521552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
522552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
523552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
524552f7358SJed Brown   PetscScalar       *ref;
525552f7358SJed Brown   PetscErrorCode    ierr;
526552f7358SJed Brown 
527552f7358SJed Brown   PetscFunctionBegin;
528552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
529552f7358SJed Brown   {
530552f7358SJed Brown     const PetscScalar x       = ref[0];
531552f7358SJed Brown     const PetscScalar y       = ref[1];
532552f7358SJed Brown     const PetscScalar z       = ref[2];
533552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
534da80777bSKarl Rupp     PetscScalar       values[9];
535da80777bSKarl Rupp 
536da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
537da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
538da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
539da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
540da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
541da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
542da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
543da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
544da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
5451aa26658SKarl Rupp 
546552f7358SJed Brown     ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
547552f7358SJed Brown   }
548552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
549552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
550552f7358SJed Brown   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
551552f7358SJed Brown   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
552552f7358SJed Brown   PetscFunctionReturn(0);
553552f7358SJed Brown }
554552f7358SJed Brown 
555552f7358SJed Brown #undef __FUNCT__
556552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Hex_Private"
557a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
558a6dfd86eSKarl Rupp {
559fafc0619SMatthew G Knepley   DM             dmCoord;
560552f7358SJed Brown   SNES           snes;
561552f7358SJed Brown   KSP            ksp;
562552f7358SJed Brown   PC             pc;
563552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
564552f7358SJed Brown   Mat            J;
565552f7358SJed Brown   PetscScalar    *a, *coords;
566552f7358SJed Brown   PetscInt       p;
567552f7358SJed Brown   PetscErrorCode ierr;
568552f7358SJed Brown 
569552f7358SJed Brown   PetscFunctionBegin;
570552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
571fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
572552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
573552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
574552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
575552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
576552f7358SJed Brown   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
577552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
578552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
579552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
580552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
581552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
582552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
5830298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
5840298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
585552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
586552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
587552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
588552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
589552f7358SJed Brown 
590552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
591552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
592552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
593*a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
594552f7358SJed Brown     PetscScalar *xi;
595cb313848SJed Brown     PetscReal    xir[3];
596552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
597552f7358SJed Brown 
598552f7358SJed Brown     /* Can make this do all points at once */
5990298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
6000adebc6cSBarry Smith     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
6010298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
6020adebc6cSBarry Smith     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
6030298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
6040298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
605552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
606552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
607552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
608552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
609552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
610552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
611552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
612cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
613cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
614cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
615552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
616552f7358SJed Brown       a[p*ctx->dof+comp] =
617cb313848SJed Brown         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
618cb313848SJed Brown         x[1*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
619cb313848SJed Brown         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
620cb313848SJed Brown         x[3*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
621cb313848SJed Brown         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
622cb313848SJed Brown         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
623cb313848SJed Brown         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
624cb313848SJed Brown         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
625552f7358SJed Brown     }
626552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
6270298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
6280298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
629552f7358SJed Brown   }
630552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
631552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
632552f7358SJed Brown 
633552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
634552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
635552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
636552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
637552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
638552f7358SJed Brown   PetscFunctionReturn(0);
639552f7358SJed Brown }
640552f7358SJed Brown 
641552f7358SJed Brown #undef __FUNCT__
642552f7358SJed Brown #define __FUNCT__ "DMInterpolationEvaluate"
643552f7358SJed Brown /*
644552f7358SJed Brown   Input Parameters:
645552f7358SJed Brown + ctx - The DMInterpolationInfo context
646552f7358SJed Brown . dm  - The DM
647552f7358SJed Brown - x   - The local vector containing the field to be interpolated
648552f7358SJed Brown 
649552f7358SJed Brown   Output Parameters:
650552f7358SJed Brown . v   - The vector containing the interpolated values
651552f7358SJed Brown */
6520adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
6530adebc6cSBarry Smith {
654552f7358SJed Brown   PetscInt       dim, coneSize, n;
655552f7358SJed Brown   PetscErrorCode ierr;
656552f7358SJed Brown 
657552f7358SJed Brown   PetscFunctionBegin;
658552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
659552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
660552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
661552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
6620adebc6cSBarry 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);
663552f7358SJed Brown   if (n) {
664552f7358SJed Brown     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
665552f7358SJed Brown     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
666552f7358SJed Brown     if (dim == 2) {
667552f7358SJed Brown       if (coneSize == 3) {
668552f7358SJed Brown         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
669552f7358SJed Brown       } else if (coneSize == 4) {
670552f7358SJed Brown         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
6710adebc6cSBarry Smith       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
672552f7358SJed Brown     } else if (dim == 3) {
673552f7358SJed Brown       if (coneSize == 4) {
674552f7358SJed Brown         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
675552f7358SJed Brown       } else {
676552f7358SJed Brown         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
677552f7358SJed Brown       }
6780adebc6cSBarry Smith     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
679552f7358SJed Brown   }
680552f7358SJed Brown   PetscFunctionReturn(0);
681552f7358SJed Brown }
682552f7358SJed Brown 
683552f7358SJed Brown #undef __FUNCT__
684552f7358SJed Brown #define __FUNCT__ "DMInterpolationDestroy"
6850adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
6860adebc6cSBarry Smith {
687552f7358SJed Brown   PetscErrorCode ierr;
688552f7358SJed Brown 
689552f7358SJed Brown   PetscFunctionBegin;
690552f7358SJed Brown   PetscValidPointer(ctx, 2);
691552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
692552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
693552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
694552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
6950298fd71SBarry Smith   *ctx = NULL;
696552f7358SJed Brown   PetscFunctionReturn(0);
697552f7358SJed Brown }
698