xref: /petsc/src/snes/utils/dmplexsnes.c (revision 3e1910f1ab6113d8365e15c6b8c907ccce7ce4ea)
1 #include <petscdmplex.h> /*I "petscdmplex.h" I*/
2 #include <petscsnes.h>      /*I "petscsnes.h" I*/
3 #include <petsc-private/petscimpl.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMInterpolationCreate"
7 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
8 {
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   PetscValidPointer(ctx, 2);
13   ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr);
14 
15   (*ctx)->comm   = comm;
16   (*ctx)->dim    = -1;
17   (*ctx)->nInput = 0;
18   (*ctx)->points = NULL;
19   (*ctx)->cells  = NULL;
20   (*ctx)->n      = -1;
21   (*ctx)->coords = NULL;
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNCT__
26 #define __FUNCT__ "DMInterpolationSetDim"
27 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
28 {
29   PetscFunctionBegin;
30   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);
31   ctx->dim = dim;
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "DMInterpolationGetDim"
37 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
38 {
39   PetscFunctionBegin;
40   PetscValidIntPointer(dim, 2);
41   *dim = ctx->dim;
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "DMInterpolationSetDof"
47 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
48 {
49   PetscFunctionBegin;
50   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);
51   ctx->dof = dof;
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "DMInterpolationGetDof"
57 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
58 {
59   PetscFunctionBegin;
60   PetscValidIntPointer(dof, 2);
61   *dof = ctx->dof;
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "DMInterpolationAddPoints"
67 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
68 {
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
73   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
74   ctx->nInput = n;
75 
76   ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr);
77   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMInterpolationSetUp"
83 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
84 {
85   MPI_Comm       comm = ctx->comm;
86   PetscScalar    *a;
87   PetscInt       p, q, i;
88   PetscMPIInt    rank, size;
89   PetscErrorCode ierr;
90   Vec            pointVec;
91   IS             cellIS;
92   PetscLayout    layout;
93   PetscReal      *globalPoints;
94   PetscScalar    *globalPointsScalar;
95   const PetscInt *ranges;
96   PetscMPIInt    *counts, *displs;
97   const PetscInt *foundCells;
98   PetscMPIInt    *foundProcs, *globalProcs;
99   PetscInt       n, N;
100 
101   PetscFunctionBegin;
102   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
103   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
104   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
105   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
106   /* Locate points */
107   n = ctx->nInput;
108   if (!redundantPoints) {
109     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
110     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
111     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
112     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
113     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
114     /* Communicate all points to all processes */
115     ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr);
116     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
117     for (p = 0; p < size; ++p) {
118       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
119       displs[p] = ranges[p]*ctx->dim;
120     }
121     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
122   } else {
123     N = n;
124 
125     globalPoints = ctx->points;
126   }
127 #if 0
128   ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
129   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
130 #else
131 #if defined(PETSC_USE_COMPLEX)
132   ierr = PetscMalloc(N*sizeof(PetscScalar),&globalPointsScalar);CHKERRQ(ierr);
133   for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i];
134 #else
135   globalPointsScalar = globalPoints;
136 #endif
137   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
138   ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
139   ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr);
140   ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr);
141 #endif
142   for (p = 0; p < N; ++p) {
143     if (foundCells[p] >= 0) foundProcs[p] = rank;
144     else foundProcs[p] = size;
145   }
146   /* Let the lowest rank process own each point */
147   ierr   = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
148   ctx->n = 0;
149   for (p = 0; p < N; ++p) {
150     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);
151     else if (globalProcs[p] == rank) ctx->n++;
152   }
153   /* Create coordinates vector and array of owned cells */
154   ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr);
155   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
156   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
157   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
158   ierr = VecSetFromOptions(ctx->coords);CHKERRQ(ierr);
159   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
160   for (p = 0, q = 0, i = 0; p < N; ++p) {
161     if (globalProcs[p] == rank) {
162       PetscInt d;
163 
164       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
165       ctx->cells[q++] = foundCells[p];
166     }
167   }
168   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
169 #if 0
170   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
171 #else
172   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
173   ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr);
174   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
175   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
176 #endif
177   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
178   if (!redundantPoints) {
179     ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);
180     ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
181   }
182   PetscFunctionReturn(0);
183 }
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "DMInterpolationGetCoordinates"
187 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
188 {
189   PetscFunctionBegin;
190   PetscValidPointer(coordinates, 2);
191   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
192   *coordinates = ctx->coords;
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "DMInterpolationGetVector"
198 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
199 {
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   PetscValidPointer(v, 2);
204   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
205   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
206   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
207   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
208   ierr = VecSetFromOptions(*v);CHKERRQ(ierr);
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "DMInterpolationRestoreVector"
214 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
215 {
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   PetscValidPointer(v, 2);
220   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
221   ierr = VecDestroy(v);CHKERRQ(ierr);
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "DMInterpolate_Simplex_Private"
227 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
228 {
229   PetscReal      *v0, *J, *invJ, detJ;
230   PetscScalar    *a, *coords;
231   PetscInt       p;
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235   ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr);
236   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
237   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
238   for (p = 0; p < ctx->n; ++p) {
239     PetscInt     c = ctx->cells[p];
240     PetscScalar *x;
241     PetscReal    xi[4];
242     PetscInt     d, f, comp;
243 
244     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
245     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
246     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
247     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
248 
249     for (d = 0; d < ctx->dim; ++d) {
250       xi[d] = 0.0;
251       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
252       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];
253     }
254     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
255   }
256   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
257   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
258   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "QuadMap_Private"
264 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
265 {
266   const PetscScalar *vertices = (const PetscScalar*) ctx;
267   const PetscScalar x0        = vertices[0];
268   const PetscScalar y0        = vertices[1];
269   const PetscScalar x1        = vertices[2];
270   const PetscScalar y1        = vertices[3];
271   const PetscScalar x2        = vertices[4];
272   const PetscScalar y2        = vertices[5];
273   const PetscScalar x3        = vertices[6];
274   const PetscScalar y3        = vertices[7];
275   const PetscScalar f_1       = x1 - x0;
276   const PetscScalar g_1       = y1 - y0;
277   const PetscScalar f_3       = x3 - x0;
278   const PetscScalar g_3       = y3 - y0;
279   const PetscScalar f_01      = x2 - x1 - x3 + x0;
280   const PetscScalar g_01      = y2 - y1 - y3 + y0;
281   PetscScalar       *ref, *real;
282   PetscErrorCode    ierr;
283 
284   PetscFunctionBegin;
285   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
286   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
287   {
288     const PetscScalar p0 = ref[0];
289     const PetscScalar p1 = ref[1];
290 
291     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
292     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
293   }
294   ierr = PetscLogFlops(28);CHKERRQ(ierr);
295   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
296   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "QuadJacobian_Private"
302 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
303 {
304   const PetscScalar *vertices = (const PetscScalar*) ctx;
305   const PetscScalar x0        = vertices[0];
306   const PetscScalar y0        = vertices[1];
307   const PetscScalar x1        = vertices[2];
308   const PetscScalar y1        = vertices[3];
309   const PetscScalar x2        = vertices[4];
310   const PetscScalar y2        = vertices[5];
311   const PetscScalar x3        = vertices[6];
312   const PetscScalar y3        = vertices[7];
313   const PetscScalar f_01      = x2 - x1 - x3 + x0;
314   const PetscScalar g_01      = y2 - y1 - y3 + y0;
315   PetscScalar       *ref;
316   PetscErrorCode    ierr;
317 
318   PetscFunctionBegin;
319   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
320   {
321     const PetscScalar x       = ref[0];
322     const PetscScalar y       = ref[1];
323     const PetscInt    rows[2] = {0, 1};
324     PetscScalar       values[4];
325 
326     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
327     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
328     ierr      = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
329   }
330   ierr = PetscLogFlops(30);CHKERRQ(ierr);
331   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
332   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334   PetscFunctionReturn(0);
335 }
336 
337 #undef __FUNCT__
338 #define __FUNCT__ "DMInterpolate_Quad_Private"
339 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
340 {
341   DM             dmCoord;
342   SNES           snes;
343   KSP            ksp;
344   PC             pc;
345   Vec            coordsLocal, r, ref, real;
346   Mat            J;
347   PetscScalar    *a, *coords;
348   PetscInt       p;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
353   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
354   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
355   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
356   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
357   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
358   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
359   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
360   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
361   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
362   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
363   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
364   ierr = MatSetUp(J);CHKERRQ(ierr);
365   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
366   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
367   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
368   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
369   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
370   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
371 
372   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
373   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
374   for (p = 0; p < ctx->n; ++p) {
375     PetscScalar *x, *vertices;
376     PetscScalar *xi;
377     PetscReal    xir[2];
378     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
379 
380     /* Can make this do all points at once */
381     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
382     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
383     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
384     if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);
385     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
386     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
387     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
388     xi[0]  = coords[p*ctx->dim+0];
389     xi[1]  = coords[p*ctx->dim+1];
390     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
391     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
392     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
393     xir[0] = PetscRealPart(xi[0]);
394     xir[1] = PetscRealPart(xi[1]);
395     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];
396 
397     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
398     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
399     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
400   }
401   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
402   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
403 
404   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
405   ierr = VecDestroy(&r);CHKERRQ(ierr);
406   ierr = VecDestroy(&ref);CHKERRQ(ierr);
407   ierr = VecDestroy(&real);CHKERRQ(ierr);
408   ierr = MatDestroy(&J);CHKERRQ(ierr);
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "HexMap_Private"
414 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
415 {
416   const PetscScalar *vertices = (const PetscScalar*) ctx;
417   const PetscScalar x0        = vertices[0];
418   const PetscScalar y0        = vertices[1];
419   const PetscScalar z0        = vertices[2];
420   const PetscScalar x1        = vertices[3];
421   const PetscScalar y1        = vertices[4];
422   const PetscScalar z1        = vertices[5];
423   const PetscScalar x2        = vertices[6];
424   const PetscScalar y2        = vertices[7];
425   const PetscScalar z2        = vertices[8];
426   const PetscScalar x3        = vertices[9];
427   const PetscScalar y3        = vertices[10];
428   const PetscScalar z3        = vertices[11];
429   const PetscScalar x4        = vertices[12];
430   const PetscScalar y4        = vertices[13];
431   const PetscScalar z4        = vertices[14];
432   const PetscScalar x5        = vertices[15];
433   const PetscScalar y5        = vertices[16];
434   const PetscScalar z5        = vertices[17];
435   const PetscScalar x6        = vertices[18];
436   const PetscScalar y6        = vertices[19];
437   const PetscScalar z6        = vertices[20];
438   const PetscScalar x7        = vertices[21];
439   const PetscScalar y7        = vertices[22];
440   const PetscScalar z7        = vertices[23];
441   const PetscScalar f_1       = x1 - x0;
442   const PetscScalar g_1       = y1 - y0;
443   const PetscScalar h_1       = z1 - z0;
444   const PetscScalar f_3       = x3 - x0;
445   const PetscScalar g_3       = y3 - y0;
446   const PetscScalar h_3       = z3 - z0;
447   const PetscScalar f_4       = x4 - x0;
448   const PetscScalar g_4       = y4 - y0;
449   const PetscScalar h_4       = z4 - z0;
450   const PetscScalar f_01      = x2 - x1 - x3 + x0;
451   const PetscScalar g_01      = y2 - y1 - y3 + y0;
452   const PetscScalar h_01      = z2 - z1 - z3 + z0;
453   const PetscScalar f_12      = x7 - x3 - x4 + x0;
454   const PetscScalar g_12      = y7 - y3 - y4 + y0;
455   const PetscScalar h_12      = z7 - z3 - z4 + z0;
456   const PetscScalar f_02      = x5 - x1 - x4 + x0;
457   const PetscScalar g_02      = y5 - y1 - y4 + y0;
458   const PetscScalar h_02      = z5 - z1 - z4 + z0;
459   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
460   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
461   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
462   PetscScalar       *ref, *real;
463   PetscErrorCode    ierr;
464 
465   PetscFunctionBegin;
466   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
467   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
468   {
469     const PetscScalar p0 = ref[0];
470     const PetscScalar p1 = ref[1];
471     const PetscScalar p2 = ref[2];
472 
473     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;
474     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;
475     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;
476   }
477   ierr = PetscLogFlops(114);CHKERRQ(ierr);
478   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
479   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "HexJacobian_Private"
485 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
486 {
487   const PetscScalar *vertices = (const PetscScalar*) ctx;
488   const PetscScalar x0        = vertices[0];
489   const PetscScalar y0        = vertices[1];
490   const PetscScalar z0        = vertices[2];
491   const PetscScalar x1        = vertices[3];
492   const PetscScalar y1        = vertices[4];
493   const PetscScalar z1        = vertices[5];
494   const PetscScalar x2        = vertices[6];
495   const PetscScalar y2        = vertices[7];
496   const PetscScalar z2        = vertices[8];
497   const PetscScalar x3        = vertices[9];
498   const PetscScalar y3        = vertices[10];
499   const PetscScalar z3        = vertices[11];
500   const PetscScalar x4        = vertices[12];
501   const PetscScalar y4        = vertices[13];
502   const PetscScalar z4        = vertices[14];
503   const PetscScalar x5        = vertices[15];
504   const PetscScalar y5        = vertices[16];
505   const PetscScalar z5        = vertices[17];
506   const PetscScalar x6        = vertices[18];
507   const PetscScalar y6        = vertices[19];
508   const PetscScalar z6        = vertices[20];
509   const PetscScalar x7        = vertices[21];
510   const PetscScalar y7        = vertices[22];
511   const PetscScalar z7        = vertices[23];
512   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
513   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
514   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
515   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
516   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
517   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
518   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
519   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
520   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
521   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
522   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
523   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
524   PetscScalar       *ref;
525   PetscErrorCode    ierr;
526 
527   PetscFunctionBegin;
528   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
529   {
530     const PetscScalar x       = ref[0];
531     const PetscScalar y       = ref[1];
532     const PetscScalar z       = ref[2];
533     const PetscInt    rows[3] = {0, 1, 2};
534     PetscScalar       values[9];
535 
536     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
537     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
538     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
539     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
540     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
541     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
542     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
543     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
544     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
545 
546     ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
547   }
548   ierr = PetscLogFlops(152);CHKERRQ(ierr);
549   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
550   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
551   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
552   PetscFunctionReturn(0);
553 }
554 
555 #undef __FUNCT__
556 #define __FUNCT__ "DMInterpolate_Hex_Private"
557 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
558 {
559   DM             dmCoord;
560   SNES           snes;
561   KSP            ksp;
562   PC             pc;
563   Vec            coordsLocal, r, ref, real;
564   Mat            J;
565   PetscScalar    *a, *coords;
566   PetscInt       p;
567   PetscErrorCode ierr;
568 
569   PetscFunctionBegin;
570   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
571   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
572   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
573   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
574   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
575   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
576   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
577   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
578   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
579   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
580   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
581   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
582   ierr = MatSetUp(J);CHKERRQ(ierr);
583   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
584   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
585   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
586   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
587   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
588   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
589 
590   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
591   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
592   for (p = 0; p < ctx->n; ++p) {
593     PetscScalar *x, *vertices;
594     PetscScalar *xi;
595     PetscReal    xir[3];
596     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
597 
598     /* Can make this do all points at once */
599     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
600     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
601     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
602     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
603     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
604     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
605     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
606     xi[0]  = coords[p*ctx->dim+0];
607     xi[1]  = coords[p*ctx->dim+1];
608     xi[2]  = coords[p*ctx->dim+2];
609     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
610     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
611     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
612     xir[0] = PetscRealPart(xi[0]);
613     xir[1] = PetscRealPart(xi[1]);
614     xir[2] = PetscRealPart(xi[2]);
615     for (comp = 0; comp < ctx->dof; ++comp) {
616       a[p*ctx->dof+comp] =
617         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
618         x[1*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
619         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
620         x[3*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
621         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
622         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
623         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
624         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
625     }
626     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
627     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
628     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
629   }
630   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
631   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
632 
633   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
634   ierr = VecDestroy(&r);CHKERRQ(ierr);
635   ierr = VecDestroy(&ref);CHKERRQ(ierr);
636   ierr = VecDestroy(&real);CHKERRQ(ierr);
637   ierr = MatDestroy(&J);CHKERRQ(ierr);
638   PetscFunctionReturn(0);
639 }
640 
641 #undef __FUNCT__
642 #define __FUNCT__ "DMInterpolationEvaluate"
643 /*
644   Input Parameters:
645 + ctx - The DMInterpolationInfo context
646 . dm  - The DM
647 - x   - The local vector containing the field to be interpolated
648 
649   Output Parameters:
650 . v   - The vector containing the interpolated values
651 */
652 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
653 {
654   PetscInt       dim, coneSize, n;
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
659   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
660   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
661   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
662   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);
663   if (n) {
664     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
665     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
666     if (dim == 2) {
667       if (coneSize == 3) {
668         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
669       } else if (coneSize == 4) {
670         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
671       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
672     } else if (dim == 3) {
673       if (coneSize == 4) {
674         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
675       } else {
676         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
677       }
678     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
679   }
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "DMInterpolationDestroy"
685 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
686 {
687   PetscErrorCode ierr;
688 
689   PetscFunctionBegin;
690   PetscValidPointer(ctx, 2);
691   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
692   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
693   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
694   ierr = PetscFree(*ctx);CHKERRQ(ierr);
695   *ctx = NULL;
696   PetscFunctionReturn(0);
697 }
698