xref: /petsc/src/snes/utils/dmplexsnes.c (revision 88c71bb6a56f75c434ed97baceb1e249d2929ce7)
1552f7358SJed Brown #include <petscdmplex.h> /*I "petscdmplex.h" I*/
2552f7358SJed Brown #include <petscsnes.h>      /*I "petscsnes.h" I*/
3552f7358SJed Brown 
4552f7358SJed Brown #undef __FUNCT__
5552f7358SJed Brown #define __FUNCT__ "DMInterpolationCreate"
6552f7358SJed Brown PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) {
7552f7358SJed Brown   PetscErrorCode ierr;
8552f7358SJed Brown 
9552f7358SJed Brown   PetscFunctionBegin;
10552f7358SJed Brown   PetscValidPointer(ctx, 2);
11552f7358SJed Brown   ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr);
12552f7358SJed Brown   (*ctx)->comm   = comm;
13552f7358SJed Brown   (*ctx)->dim    = -1;
14552f7358SJed Brown   (*ctx)->nInput = 0;
15552f7358SJed Brown   (*ctx)->points = PETSC_NULL;
16552f7358SJed Brown   (*ctx)->cells  = PETSC_NULL;
17552f7358SJed Brown   (*ctx)->n      = -1;
18552f7358SJed Brown   (*ctx)->coords = PETSC_NULL;
19552f7358SJed Brown   PetscFunctionReturn(0);
20552f7358SJed Brown }
21552f7358SJed Brown 
22552f7358SJed Brown #undef __FUNCT__
23552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDim"
24552f7358SJed Brown PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) {
25552f7358SJed Brown   PetscFunctionBegin;
26552f7358SJed Brown   if ((dim < 1) || (dim > 3)) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);}
27552f7358SJed Brown   ctx->dim = dim;
28552f7358SJed Brown   PetscFunctionReturn(0);
29552f7358SJed Brown }
30552f7358SJed Brown 
31552f7358SJed Brown #undef __FUNCT__
32552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDim"
33552f7358SJed Brown PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) {
34552f7358SJed Brown   PetscFunctionBegin;
35552f7358SJed Brown   PetscValidIntPointer(dim, 2);
36552f7358SJed Brown   *dim = ctx->dim;
37552f7358SJed Brown   PetscFunctionReturn(0);
38552f7358SJed Brown }
39552f7358SJed Brown 
40552f7358SJed Brown #undef __FUNCT__
41552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDof"
42552f7358SJed Brown PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) {
43552f7358SJed Brown   PetscFunctionBegin;
44552f7358SJed Brown   if (dof < 1) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);}
45552f7358SJed Brown   ctx->dof = dof;
46552f7358SJed Brown   PetscFunctionReturn(0);
47552f7358SJed Brown }
48552f7358SJed Brown 
49552f7358SJed Brown #undef __FUNCT__
50552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDof"
51552f7358SJed Brown PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) {
52552f7358SJed Brown   PetscFunctionBegin;
53552f7358SJed Brown   PetscValidIntPointer(dof, 2);
54552f7358SJed Brown   *dof = ctx->dof;
55552f7358SJed Brown   PetscFunctionReturn(0);
56552f7358SJed Brown }
57552f7358SJed Brown 
58552f7358SJed Brown #undef __FUNCT__
59552f7358SJed Brown #define __FUNCT__ "DMInterpolationAddPoints"
60552f7358SJed Brown PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) {
61552f7358SJed Brown   PetscErrorCode ierr;
62552f7358SJed Brown 
63552f7358SJed Brown   PetscFunctionBegin;
64552f7358SJed Brown   if (ctx->dim < 0) {
65552f7358SJed Brown     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
66552f7358SJed Brown   }
67552f7358SJed Brown   if (ctx->points) {
68552f7358SJed Brown     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
69552f7358SJed Brown   }
70552f7358SJed Brown   ctx->nInput = n;
71552f7358SJed Brown   ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr);
72552f7358SJed Brown   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
73552f7358SJed Brown   PetscFunctionReturn(0);
74552f7358SJed Brown }
75552f7358SJed Brown 
76552f7358SJed Brown #undef __FUNCT__
77552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetUp"
78552f7358SJed Brown PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) {
79552f7358SJed Brown   MPI_Comm       comm = ctx->comm;
80552f7358SJed Brown   PetscScalar   *a;
81552f7358SJed Brown   PetscInt       p, q, i;
82552f7358SJed Brown   PetscMPIInt    rank, size;
83552f7358SJed Brown   PetscErrorCode ierr;
84552f7358SJed Brown 
85552f7358SJed Brown   PetscFunctionBegin;
86552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87552f7358SJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
88552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
89552f7358SJed Brown   if (ctx->dim < 0) {
90552f7358SJed Brown     SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
91552f7358SJed Brown   }
92552f7358SJed Brown   // Locate points
93552f7358SJed Brown   Vec             pointVec;
94552f7358SJed Brown   IS              cellIS;
95552f7358SJed Brown   PetscLayout     layout;
96552f7358SJed Brown   PetscReal      *globalPoints;
97552f7358SJed Brown   const PetscInt *ranges;
98552f7358SJed Brown   PetscMPIInt    *counts, *displs;
99552f7358SJed Brown   const PetscInt *foundCells;
100552f7358SJed Brown   PetscMPIInt    *foundProcs, *globalProcs;
101552f7358SJed Brown   PetscInt        n = ctx->nInput, N;
102552f7358SJed Brown 
103552f7358SJed Brown   if (!redundantPoints) {
104552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
105552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
106552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
107552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
108552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
109552f7358SJed Brown     /* Communicate all points to all processes */
110552f7358SJed Brown     ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr);
111552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
112552f7358SJed Brown     for (p = 0; p < size; ++p) {
113552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
114552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
115552f7358SJed Brown     }
116552f7358SJed Brown     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
117552f7358SJed Brown   } else {
118552f7358SJed Brown     N = n;
119552f7358SJed Brown     globalPoints = ctx->points;
120552f7358SJed Brown   }
121552f7358SJed Brown #if 0
122552f7358SJed Brown   ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
123552f7358SJed Brown   //foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]);
124552f7358SJed Brown #else
125552f7358SJed Brown   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N, globalPoints, &pointVec);CHKERRQ(ierr);
126552f7358SJed Brown   ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
127552f7358SJed Brown   ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr);
128552f7358SJed Brown   ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr);
129552f7358SJed Brown #endif
130552f7358SJed Brown   for (p = 0; p < N; ++p) {
131552f7358SJed Brown     if (foundCells[p] >= 0) {
132552f7358SJed Brown       foundProcs[p] = rank;
133552f7358SJed Brown     } else {
134552f7358SJed Brown       foundProcs[p] = size;
135552f7358SJed Brown     }
136552f7358SJed Brown   }
137552f7358SJed Brown   /* Let the lowest rank process own each point */
138552f7358SJed Brown   ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
139552f7358SJed Brown   ctx->n = 0;
140552f7358SJed Brown   for (p = 0; p < N; ++p) {
141552f7358SJed Brown     if (globalProcs[p] == size) {
142552f7358SJed Brown       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);
143552f7358SJed Brown     } else if (globalProcs[p] == rank) {
144552f7358SJed Brown       ctx->n++;
145552f7358SJed Brown     }
146552f7358SJed Brown   }
147552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
148552f7358SJed Brown   ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr);
149552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
150552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
151552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
152552f7358SJed Brown   ierr = VecSetFromOptions(ctx->coords);CHKERRQ(ierr);
153552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
154552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
155552f7358SJed Brown     if (globalProcs[p] == rank) {
156552f7358SJed Brown       PetscInt d;
157552f7358SJed Brown 
158552f7358SJed Brown       for (d = 0; d < ctx->dim; ++d, ++i) {
159552f7358SJed Brown         a[i] = globalPoints[p*ctx->dim+d];
160552f7358SJed Brown       }
161552f7358SJed Brown       ctx->cells[q++] = foundCells[p];
162552f7358SJed Brown     }
163552f7358SJed Brown   }
164552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
165552f7358SJed Brown #if 0
166552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
167552f7358SJed Brown #else
168552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
169552f7358SJed Brown   ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr);
170552f7358SJed Brown   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
171552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
172552f7358SJed Brown #endif
173552f7358SJed Brown   if (!redundantPoints) {
174552f7358SJed Brown     ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);
175552f7358SJed Brown     ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
176552f7358SJed Brown   }
177552f7358SJed Brown   PetscFunctionReturn(0);
178552f7358SJed Brown }
179552f7358SJed Brown 
180552f7358SJed Brown #undef __FUNCT__
181552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetCoordinates"
182552f7358SJed Brown PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) {
183552f7358SJed Brown   PetscFunctionBegin;
184552f7358SJed Brown   PetscValidPointer(coordinates, 2);
185552f7358SJed Brown   if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
186552f7358SJed Brown   *coordinates = ctx->coords;
187552f7358SJed Brown   PetscFunctionReturn(0);
188552f7358SJed Brown }
189552f7358SJed Brown 
190552f7358SJed Brown #undef __FUNCT__
191552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetVector"
192552f7358SJed Brown PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) {
193552f7358SJed Brown   PetscErrorCode ierr;
194552f7358SJed Brown 
195552f7358SJed Brown   PetscFunctionBegin;
196552f7358SJed Brown   PetscValidPointer(v, 2);
197552f7358SJed Brown   if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
198552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
199552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
200552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
201552f7358SJed Brown   ierr = VecSetFromOptions(*v);CHKERRQ(ierr);
202552f7358SJed Brown   PetscFunctionReturn(0);
203552f7358SJed Brown }
204552f7358SJed Brown 
205552f7358SJed Brown #undef __FUNCT__
206552f7358SJed Brown #define __FUNCT__ "DMInterpolationRestoreVector"
207552f7358SJed Brown PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) {
208552f7358SJed Brown   PetscErrorCode ierr;
209552f7358SJed Brown 
210552f7358SJed Brown   PetscFunctionBegin;
211552f7358SJed Brown   PetscValidPointer(v, 2);
212552f7358SJed Brown   if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
213552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
214552f7358SJed Brown   PetscFunctionReturn(0);
215552f7358SJed Brown }
216552f7358SJed Brown 
217552f7358SJed Brown #undef __FUNCT__
218552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Simplex_Private"
219552f7358SJed Brown PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
220552f7358SJed Brown   PetscReal     *v0, *J, *invJ, detJ;
221552f7358SJed Brown   PetscScalar   *a, *coords;
222552f7358SJed Brown   PetscInt       p;
223552f7358SJed Brown   PetscErrorCode ierr;
224552f7358SJed Brown 
225552f7358SJed Brown   PetscFunctionBegin;
226552f7358SJed Brown   ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr);
227552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
228552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
229552f7358SJed Brown   for(p = 0; p < ctx->n; ++p) {
230552f7358SJed Brown     PetscInt           c = ctx->cells[p];
231552f7358SJed Brown     const PetscScalar *x;
232552f7358SJed Brown     PetscReal          xi[4];
233552f7358SJed Brown     PetscInt           d, f, comp;
234552f7358SJed Brown 
235552f7358SJed Brown     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
236552f7358SJed Brown     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
237552f7358SJed Brown     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr);
238552f7358SJed Brown     for(comp = 0; comp < ctx->dof; ++comp) {
239552f7358SJed Brown       a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
240552f7358SJed Brown     }
241552f7358SJed Brown     for(d = 0; d < ctx->dim; ++d) {
242552f7358SJed Brown       xi[d] = 0.0;
243552f7358SJed Brown       for(f = 0; f < ctx->dim; ++f) {
244552f7358SJed Brown         xi[d] += invJ[d*ctx->dim+f]*0.5*(coords[p*ctx->dim+f] - v0[f]);
245552f7358SJed Brown       }
246552f7358SJed Brown       for(comp = 0; comp < ctx->dof; ++comp) {
247552f7358SJed Brown         a[p*ctx->dof+comp] += (x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
248552f7358SJed Brown       }
249552f7358SJed Brown     }
250552f7358SJed Brown     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr);
251552f7358SJed Brown   }
252552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
253552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
254552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
255552f7358SJed Brown   PetscFunctionReturn(0);
256552f7358SJed Brown }
257552f7358SJed Brown 
258552f7358SJed Brown #undef __FUNCT__
259552f7358SJed Brown #define __FUNCT__ "QuadMap_Private"
260552f7358SJed Brown PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
261552f7358SJed Brown {
262552f7358SJed Brown   const PetscScalar*vertices = (const PetscScalar *) ctx;
263552f7358SJed Brown   const PetscScalar x0   = vertices[0];
264552f7358SJed Brown   const PetscScalar y0   = vertices[1];
265552f7358SJed Brown   const PetscScalar x1   = vertices[2];
266552f7358SJed Brown   const PetscScalar y1   = vertices[3];
267552f7358SJed Brown   const PetscScalar x2   = vertices[4];
268552f7358SJed Brown   const PetscScalar y2   = vertices[5];
269552f7358SJed Brown   const PetscScalar x3   = vertices[6];
270552f7358SJed Brown   const PetscScalar y3   = vertices[7];
271552f7358SJed Brown   const PetscScalar f_1  = x1 - x0;
272552f7358SJed Brown   const PetscScalar g_1  = y1 - y0;
273552f7358SJed Brown   const PetscScalar f_3  = x3 - x0;
274552f7358SJed Brown   const PetscScalar g_3  = y3 - y0;
275552f7358SJed Brown   const PetscScalar f_01 = x2 - x1 - x3 + x0;
276552f7358SJed Brown   const PetscScalar g_01 = y2 - y1 - y3 + y0;
277552f7358SJed Brown   PetscScalar      *ref, *real;
278552f7358SJed Brown   PetscErrorCode    ierr;
279552f7358SJed Brown 
280552f7358SJed Brown   PetscFunctionBegin;
281552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
282552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
283552f7358SJed Brown   {
284552f7358SJed Brown     const PetscScalar p0 = ref[0];
285552f7358SJed Brown     const PetscScalar p1 = ref[1];
286552f7358SJed Brown 
287552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
288552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
289552f7358SJed Brown   }
290552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
291552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
292552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
293552f7358SJed Brown   PetscFunctionReturn(0);
294552f7358SJed Brown }
295552f7358SJed Brown 
296552f7358SJed Brown #undef __FUNCT__
297552f7358SJed Brown #define __FUNCT__ "QuadJacobian_Private"
298552f7358SJed Brown PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
299552f7358SJed Brown {
300552f7358SJed Brown   const PetscScalar*vertices = (const PetscScalar *) ctx;
301552f7358SJed Brown   const PetscScalar x0   = vertices[0];
302552f7358SJed Brown   const PetscScalar y0   = vertices[1];
303552f7358SJed Brown   const PetscScalar x1   = vertices[2];
304552f7358SJed Brown   const PetscScalar y1   = vertices[3];
305552f7358SJed Brown   const PetscScalar x2   = vertices[4];
306552f7358SJed Brown   const PetscScalar y2   = vertices[5];
307552f7358SJed Brown   const PetscScalar x3   = vertices[6];
308552f7358SJed Brown   const PetscScalar y3   = vertices[7];
309552f7358SJed Brown   const PetscScalar f_01 = x2 - x1 - x3 + x0;
310552f7358SJed Brown   const PetscScalar g_01 = y2 - y1 - y3 + y0;
311552f7358SJed Brown   PetscScalar      *ref;
312552f7358SJed Brown   PetscErrorCode    ierr;
313552f7358SJed Brown 
314552f7358SJed Brown   PetscFunctionBegin;
315552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
316552f7358SJed Brown   {
317552f7358SJed Brown     const PetscScalar x = ref[0];
318552f7358SJed Brown     const PetscScalar y = ref[1];
319552f7358SJed Brown     const PetscInt    rows[2]   = {0, 1};
320552f7358SJed Brown     const PetscScalar values[4] = {(x1 - x0 + f_01*y) * 0.5, (x3 - x0 + f_01*x) * 0.5,
321552f7358SJed Brown                                    (y1 - y0 + g_01*y) * 0.5, (y3 - y0 + g_01*x) * 0.5};
322552f7358SJed Brown     ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
323552f7358SJed Brown   }
324552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
325552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
326552f7358SJed Brown   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327552f7358SJed Brown   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
328552f7358SJed Brown   PetscFunctionReturn(0);
329552f7358SJed Brown }
330552f7358SJed Brown 
331552f7358SJed Brown #undef __FUNCT__
332552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Quad_Private"
333552f7358SJed Brown PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
334552f7358SJed Brown   SNES           snes;
335552f7358SJed Brown   KSP            ksp;
336552f7358SJed Brown   PC             pc;
337552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
338552f7358SJed Brown   Mat            J;
339552f7358SJed Brown   PetscScalar   *a, *coords;
340552f7358SJed Brown   PetscInt       p;
341552f7358SJed Brown   PetscErrorCode ierr;
342552f7358SJed Brown 
343552f7358SJed Brown   PetscFunctionBegin;
344552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
345552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
346552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
347552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
348552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
349552f7358SJed Brown   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
350552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
351552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
352552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
353552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
354552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
355552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
356552f7358SJed Brown   ierr = SNESSetFunction(snes, r, QuadMap_Private, PETSC_NULL);CHKERRQ(ierr);
357552f7358SJed Brown   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, PETSC_NULL);CHKERRQ(ierr);
358552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
359552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
360552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
361552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
362552f7358SJed Brown 
363552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
364552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
365552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
366552f7358SJed Brown     const PetscScalar *x, *vertices;
367552f7358SJed Brown     PetscScalar       *xi;
368552f7358SJed Brown     PetscInt           c = ctx->cells[p], comp, coordSize, xSize;
369552f7358SJed Brown 
370552f7358SJed Brown     /* Can make this do all points at once */
371552f7358SJed Brown     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
372552f7358SJed Brown     if (4*2 != coordSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);}
373552f7358SJed Brown     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
374552f7358SJed Brown     if (4*ctx->dof != xSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);}
375*88c71bb6SMatthew G Knepley     ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
376*88c71bb6SMatthew G Knepley     ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
377552f7358SJed Brown     ierr = VecGetArray(real, &xi);CHKERRQ(ierr);
378552f7358SJed Brown     xi[0] = coords[p*ctx->dim+0];
379552f7358SJed Brown     xi[1] = coords[p*ctx->dim+1];
380552f7358SJed Brown     ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr);
381552f7358SJed Brown     ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr);
382552f7358SJed Brown     ierr = VecGetArray(ref, &xi);CHKERRQ(ierr);
383552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
384552f7358SJed Brown       a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xi[0])*(1 - xi[1]) + x[1*ctx->dof+comp]*xi[0]*(1 - xi[1]) + x[2*ctx->dof+comp]*xi[0]*xi[1] + x[3*ctx->dof+comp]*(1 - xi[0])*xi[1];
385552f7358SJed Brown     }
386552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
387552f7358SJed Brown     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
388552f7358SJed Brown   }
389552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
390552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
391552f7358SJed Brown 
392552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
393552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
394552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
395552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
396552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
397552f7358SJed Brown   PetscFunctionReturn(0);
398552f7358SJed Brown }
399552f7358SJed Brown 
400552f7358SJed Brown #undef __FUNCT__
401552f7358SJed Brown #define __FUNCT__ "HexMap_Private"
402552f7358SJed Brown PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
403552f7358SJed Brown {
404552f7358SJed Brown   const PetscScalar*vertices = (const PetscScalar *) ctx;
405552f7358SJed Brown   const PetscScalar x0 = vertices[0];
406552f7358SJed Brown   const PetscScalar y0 = vertices[1];
407552f7358SJed Brown   const PetscScalar z0 = vertices[2];
408552f7358SJed Brown   const PetscScalar x1 = vertices[3];
409552f7358SJed Brown   const PetscScalar y1 = vertices[4];
410552f7358SJed Brown   const PetscScalar z1 = vertices[5];
411552f7358SJed Brown   const PetscScalar x2 = vertices[6];
412552f7358SJed Brown   const PetscScalar y2 = vertices[7];
413552f7358SJed Brown   const PetscScalar z2 = vertices[8];
414552f7358SJed Brown   const PetscScalar x3 = vertices[9];
415552f7358SJed Brown   const PetscScalar y3 = vertices[10];
416552f7358SJed Brown   const PetscScalar z3 = vertices[11];
417552f7358SJed Brown   const PetscScalar x4 = vertices[12];
418552f7358SJed Brown   const PetscScalar y4 = vertices[13];
419552f7358SJed Brown   const PetscScalar z4 = vertices[14];
420552f7358SJed Brown   const PetscScalar x5 = vertices[15];
421552f7358SJed Brown   const PetscScalar y5 = vertices[16];
422552f7358SJed Brown   const PetscScalar z5 = vertices[17];
423552f7358SJed Brown   const PetscScalar x6 = vertices[18];
424552f7358SJed Brown   const PetscScalar y6 = vertices[19];
425552f7358SJed Brown   const PetscScalar z6 = vertices[20];
426552f7358SJed Brown   const PetscScalar x7 = vertices[21];
427552f7358SJed Brown   const PetscScalar y7 = vertices[22];
428552f7358SJed Brown   const PetscScalar z7 = vertices[23];
429552f7358SJed Brown   const PetscScalar f_1 = x1 - x0;
430552f7358SJed Brown   const PetscScalar g_1 = y1 - y0;
431552f7358SJed Brown   const PetscScalar h_1 = z1 - z0;
432552f7358SJed Brown   const PetscScalar f_3 = x3 - x0;
433552f7358SJed Brown   const PetscScalar g_3 = y3 - y0;
434552f7358SJed Brown   const PetscScalar h_3 = z3 - z0;
435552f7358SJed Brown   const PetscScalar f_4 = x4 - x0;
436552f7358SJed Brown   const PetscScalar g_4 = y4 - y0;
437552f7358SJed Brown   const PetscScalar h_4 = z4 - z0;
438552f7358SJed Brown   const PetscScalar f_01 = x2 - x1 - x3 + x0;
439552f7358SJed Brown   const PetscScalar g_01 = y2 - y1 - y3 + y0;
440552f7358SJed Brown   const PetscScalar h_01 = z2 - z1 - z3 + z0;
441552f7358SJed Brown   const PetscScalar f_12 = x7 - x3 - x4 + x0;
442552f7358SJed Brown   const PetscScalar g_12 = y7 - y3 - y4 + y0;
443552f7358SJed Brown   const PetscScalar h_12 = z7 - z3 - z4 + z0;
444552f7358SJed Brown   const PetscScalar f_02 = x5 - x1 - x4 + x0;
445552f7358SJed Brown   const PetscScalar g_02 = y5 - y1 - y4 + y0;
446552f7358SJed Brown   const PetscScalar h_02 = z5 - z1 - z4 + z0;
447552f7358SJed Brown   const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
448552f7358SJed Brown   const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
449552f7358SJed Brown   const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
450552f7358SJed Brown   PetscScalar      *ref, *real;
451552f7358SJed Brown   PetscErrorCode    ierr;
452552f7358SJed Brown 
453552f7358SJed Brown   PetscFunctionBegin;
454552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
455552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
456552f7358SJed Brown   {
457552f7358SJed Brown     const PetscScalar p0 = ref[0];
458552f7358SJed Brown     const PetscScalar p1 = ref[1];
459552f7358SJed Brown     const PetscScalar p2 = ref[2];
460552f7358SJed Brown 
461552f7358SJed 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;
462552f7358SJed 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;
463552f7358SJed 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;
464552f7358SJed Brown   }
465552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
466552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
467552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
468552f7358SJed Brown   PetscFunctionReturn(0);
469552f7358SJed Brown }
470552f7358SJed Brown 
471552f7358SJed Brown #undef __FUNCT__
472552f7358SJed Brown #define __FUNCT__ "HexJacobian_Private"
473552f7358SJed Brown PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
474552f7358SJed Brown {
475552f7358SJed Brown   const PetscScalar*vertices = (const PetscScalar *) ctx;
476552f7358SJed Brown   const PetscScalar x0 = vertices[0];
477552f7358SJed Brown   const PetscScalar y0 = vertices[1];
478552f7358SJed Brown   const PetscScalar z0 = vertices[2];
479552f7358SJed Brown   const PetscScalar x1 = vertices[3];
480552f7358SJed Brown   const PetscScalar y1 = vertices[4];
481552f7358SJed Brown   const PetscScalar z1 = vertices[5];
482552f7358SJed Brown   const PetscScalar x2 = vertices[6];
483552f7358SJed Brown   const PetscScalar y2 = vertices[7];
484552f7358SJed Brown   const PetscScalar z2 = vertices[8];
485552f7358SJed Brown   const PetscScalar x3 = vertices[9];
486552f7358SJed Brown   const PetscScalar y3 = vertices[10];
487552f7358SJed Brown   const PetscScalar z3 = vertices[11];
488552f7358SJed Brown   const PetscScalar x4 = vertices[12];
489552f7358SJed Brown   const PetscScalar y4 = vertices[13];
490552f7358SJed Brown   const PetscScalar z4 = vertices[14];
491552f7358SJed Brown   const PetscScalar x5 = vertices[15];
492552f7358SJed Brown   const PetscScalar y5 = vertices[16];
493552f7358SJed Brown   const PetscScalar z5 = vertices[17];
494552f7358SJed Brown   const PetscScalar x6 = vertices[18];
495552f7358SJed Brown   const PetscScalar y6 = vertices[19];
496552f7358SJed Brown   const PetscScalar z6 = vertices[20];
497552f7358SJed Brown   const PetscScalar x7 = vertices[21];
498552f7358SJed Brown   const PetscScalar y7 = vertices[22];
499552f7358SJed Brown   const PetscScalar z7 = vertices[23];
500552f7358SJed Brown   const PetscScalar f_xy = x2 - x1 - x3 + x0;
501552f7358SJed Brown   const PetscScalar g_xy = y2 - y1 - y3 + y0;
502552f7358SJed Brown   const PetscScalar h_xy = z2 - z1 - z3 + z0;
503552f7358SJed Brown   const PetscScalar f_yz = x7 - x3 - x4 + x0;
504552f7358SJed Brown   const PetscScalar g_yz = y7 - y3 - y4 + y0;
505552f7358SJed Brown   const PetscScalar h_yz = z7 - z3 - z4 + z0;
506552f7358SJed Brown   const PetscScalar f_xz = x5 - x1 - x4 + x0;
507552f7358SJed Brown   const PetscScalar g_xz = y5 - y1 - y4 + y0;
508552f7358SJed Brown   const PetscScalar h_xz = z5 - z1 - z4 + z0;
509552f7358SJed Brown   const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
510552f7358SJed Brown   const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
511552f7358SJed Brown   const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
512552f7358SJed Brown   PetscScalar      *ref;
513552f7358SJed Brown   PetscErrorCode    ierr;
514552f7358SJed Brown 
515552f7358SJed Brown   PetscFunctionBegin;
516552f7358SJed Brown   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
517552f7358SJed Brown   {
518552f7358SJed Brown     const PetscScalar x = ref[0];
519552f7358SJed Brown     const PetscScalar y = ref[1];
520552f7358SJed Brown     const PetscScalar z = ref[2];
521552f7358SJed Brown     const PetscInt    rows[3]   = {0, 1, 2};
522552f7358SJed Brown     const PetscScalar values[9] = {
523552f7358SJed Brown       (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0,
524552f7358SJed Brown       (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0,
525552f7358SJed Brown       (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0,
526552f7358SJed Brown       (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0,
527552f7358SJed Brown       (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0,
528552f7358SJed Brown       (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0,
529552f7358SJed Brown       (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0,
530552f7358SJed Brown       (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0,
531552f7358SJed Brown       (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0};
532552f7358SJed Brown     ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
533552f7358SJed Brown   }
534552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
535552f7358SJed Brown   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
536552f7358SJed Brown   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
537552f7358SJed Brown   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
538552f7358SJed Brown   PetscFunctionReturn(0);
539552f7358SJed Brown }
540552f7358SJed Brown 
541552f7358SJed Brown #undef __FUNCT__
542552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Hex_Private"
543552f7358SJed Brown PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
544552f7358SJed Brown   SNES           snes;
545552f7358SJed Brown   KSP            ksp;
546552f7358SJed Brown   PC             pc;
547552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
548552f7358SJed Brown   Mat            J;
549552f7358SJed Brown   PetscScalar   *a, *coords;
550552f7358SJed Brown   PetscInt       p;
551552f7358SJed Brown   PetscErrorCode ierr;
552552f7358SJed Brown 
553552f7358SJed Brown   PetscFunctionBegin;
554552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
555552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
556552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
557552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
558552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
559552f7358SJed Brown   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
560552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
561552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
562552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
563552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
564552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
565552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
566552f7358SJed Brown   ierr = SNESSetFunction(snes, r, HexMap_Private, PETSC_NULL);CHKERRQ(ierr);
567552f7358SJed Brown   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, PETSC_NULL);CHKERRQ(ierr);
568552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
569552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
570552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
571552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
572552f7358SJed Brown 
573552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
574552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
575552f7358SJed Brown   for(p = 0; p < ctx->n; ++p) {
576552f7358SJed Brown     const PetscScalar *x, *vertices;
577552f7358SJed Brown     PetscScalar       *xi;
578552f7358SJed Brown     PetscInt           c = ctx->cells[p], comp, coordSize, xSize;
579552f7358SJed Brown 
580552f7358SJed Brown     /* Can make this do all points at once */
581552f7358SJed Brown     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
582552f7358SJed Brown     if (8*3 != coordSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);}
583552f7358SJed Brown     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
584552f7358SJed Brown     if (8*ctx->dof != xSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);}
585*88c71bb6SMatthew G Knepley     ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
586*88c71bb6SMatthew G Knepley     ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
587552f7358SJed Brown     ierr = VecGetArray(real, &xi);CHKERRQ(ierr);
588552f7358SJed Brown     xi[0] = coords[p*ctx->dim+0];
589552f7358SJed Brown     xi[1] = coords[p*ctx->dim+1];
590552f7358SJed Brown     xi[2] = coords[p*ctx->dim+2];
591552f7358SJed Brown     ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr);
592552f7358SJed Brown     ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr);
593552f7358SJed Brown     ierr = VecGetArray(ref, &xi);CHKERRQ(ierr);
594552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
595552f7358SJed Brown       a[p*ctx->dof+comp] =
596552f7358SJed Brown         x[0*ctx->dof+comp]*(1-xi[0])*(1-xi[1])*(1-xi[2]) +
597552f7358SJed Brown         x[1*ctx->dof+comp]*    xi[0]*(1-xi[1])*(1-xi[2]) +
598552f7358SJed Brown         x[2*ctx->dof+comp]*    xi[0]*    xi[1]*(1-xi[2]) +
599552f7358SJed Brown         x[3*ctx->dof+comp]*(1-xi[0])*    xi[1]*(1-xi[2]) +
600552f7358SJed Brown         x[4*ctx->dof+comp]*(1-xi[0])*(1-xi[1])*   xi[2] +
601552f7358SJed Brown         x[5*ctx->dof+comp]*    xi[0]*(1-xi[1])*   xi[2] +
602552f7358SJed Brown         x[6*ctx->dof+comp]*    xi[0]*    xi[1]*   xi[2] +
603552f7358SJed Brown         x[7*ctx->dof+comp]*(1-xi[0])*    xi[1]*   xi[2];
604552f7358SJed Brown     }
605552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
606552f7358SJed Brown     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
607552f7358SJed Brown     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
608552f7358SJed Brown   }
609552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
610552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
611552f7358SJed Brown 
612552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
613552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
614552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
615552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
616552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
617552f7358SJed Brown   PetscFunctionReturn(0);
618552f7358SJed Brown }
619552f7358SJed Brown 
620552f7358SJed Brown #undef __FUNCT__
621552f7358SJed Brown #define __FUNCT__ "DMInterpolationEvaluate"
622552f7358SJed Brown /*
623552f7358SJed Brown   Input Parameters:
624552f7358SJed Brown + ctx - The DMInterpolationInfo context
625552f7358SJed Brown . dm  - The DM
626552f7358SJed Brown - x   - The local vector containing the field to be interpolated
627552f7358SJed Brown 
628552f7358SJed Brown   Output Parameters:
629552f7358SJed Brown . v   - The vector containing the interpolated values
630552f7358SJed Brown */
631552f7358SJed Brown PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) {
632552f7358SJed Brown   PetscInt       dim, coneSize, n;
633552f7358SJed Brown   PetscErrorCode ierr;
634552f7358SJed Brown 
635552f7358SJed Brown   PetscFunctionBegin;
636552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
637552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
638552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
639552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
640552f7358SJed Brown   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);}
641552f7358SJed Brown   if (n) {
642552f7358SJed Brown     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
643552f7358SJed Brown     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
644552f7358SJed Brown     if (dim == 2) {
645552f7358SJed Brown       if (coneSize == 3) {
646552f7358SJed Brown         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
647552f7358SJed Brown       } else if (coneSize == 4) {
648552f7358SJed Brown         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
649552f7358SJed Brown       } else {
650552f7358SJed Brown         SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
651552f7358SJed Brown       }
652552f7358SJed Brown     } else if (dim == 3) {
653552f7358SJed Brown       if (coneSize == 4) {
654552f7358SJed Brown         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
655552f7358SJed Brown       } else {
656552f7358SJed Brown         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
657552f7358SJed Brown       }
658552f7358SJed Brown     } else {
659552f7358SJed Brown       SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
660552f7358SJed Brown     }
661552f7358SJed Brown   }
662552f7358SJed Brown   PetscFunctionReturn(0);
663552f7358SJed Brown }
664552f7358SJed Brown 
665552f7358SJed Brown #undef __FUNCT__
666552f7358SJed Brown #define __FUNCT__ "DMInterpolationDestroy"
667552f7358SJed Brown PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) {
668552f7358SJed Brown   PetscErrorCode ierr;
669552f7358SJed Brown 
670552f7358SJed Brown   PetscFunctionBegin;
671552f7358SJed Brown   PetscValidPointer(ctx, 2);
672552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
673552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
674552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
675552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
676552f7358SJed Brown   *ctx = PETSC_NULL;
677552f7358SJed Brown   PetscFunctionReturn(0);
678552f7358SJed Brown }
679