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