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