xref: /petsc/src/snes/utils/dmplexsnes.c (revision 3bf036e263fd78807e2931ff42d9ddcd8aae3fd4)
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   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   PetscValidPointer(ctx, 2);
11   ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr);
12   (*ctx)->comm   = comm;
13   (*ctx)->dim    = -1;
14   (*ctx)->nInput = 0;
15   (*ctx)->points = PETSC_NULL;
16   (*ctx)->cells  = PETSC_NULL;
17   (*ctx)->n      = -1;
18   (*ctx)->coords = PETSC_NULL;
19   PetscFunctionReturn(0);
20 }
21 
22 #undef __FUNCT__
23 #define __FUNCT__ "DMInterpolationSetDim"
24 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) {
25   PetscFunctionBegin;
26   if ((dim < 1) || (dim > 3)) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);}
27   ctx->dim = dim;
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "DMInterpolationGetDim"
33 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) {
34   PetscFunctionBegin;
35   PetscValidIntPointer(dim, 2);
36   *dim = ctx->dim;
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "DMInterpolationSetDof"
42 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) {
43   PetscFunctionBegin;
44   if (dof < 1) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);}
45   ctx->dof = dof;
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "DMInterpolationGetDof"
51 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) {
52   PetscFunctionBegin;
53   PetscValidIntPointer(dof, 2);
54   *dof = ctx->dof;
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "DMInterpolationAddPoints"
60 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) {
61   PetscErrorCode ierr;
62 
63   PetscFunctionBegin;
64   if (ctx->dim < 0) {
65     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
66   }
67   if (ctx->points) {
68     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
69   }
70   ctx->nInput = n;
71   ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr);
72   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "DMInterpolationSetUp"
78 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) {
79   MPI_Comm       comm = ctx->comm;
80   PetscScalar   *a;
81   PetscInt       p, q, i;
82   PetscMPIInt    rank, size;
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
88   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
89   if (ctx->dim < 0) {
90     SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
91   }
92   // Locate points
93   Vec             pointVec;
94   IS              cellIS;
95   PetscLayout     layout;
96   PetscReal      *globalPoints;
97   PetscScalar    *globalPointsScalar;
98   const PetscInt *ranges;
99   PetscMPIInt    *counts, *displs;
100   const PetscInt *foundCells;
101   PetscMPIInt    *foundProcs, *globalProcs;
102   PetscInt        n = ctx->nInput, N;
103 
104   if (!redundantPoints) {
105     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
106     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
107     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
108     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
109     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
110     /* Communicate all points to all processes */
111     ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr);
112     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
113     for (p = 0; p < size; ++p) {
114       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
115       displs[p] = ranges[p]*ctx->dim;
116     }
117     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
118   } else {
119     N = n;
120     globalPoints = ctx->points;
121   }
122 #if 0
123   ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
124   //foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]);
125 #else
126 #if defined(PETSC_USE_COMPLEX)
127   ierr = PetscMalloc(N*sizeof(PetscScalar),&globalPointsScalar);CHKERRQ(ierr);
128   for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i];
129 #else
130   globalPointsScalar = globalPoints;
131 #endif
132   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
133   ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr);
134   ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr);
135   ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr);
136 #endif
137   for (p = 0; p < N; ++p) {
138     if (foundCells[p] >= 0) {
139       foundProcs[p] = rank;
140     } else {
141       foundProcs[p] = size;
142     }
143   }
144   /* Let the lowest rank process own each point */
145   ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPIU_INT, MPI_MIN, comm);CHKERRQ(ierr);
146   ctx->n = 0;
147   for (p = 0; p < N; ++p) {
148     if (globalProcs[p] == size) {
149       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   PetscFunctionBegin;
192   PetscValidPointer(coordinates, 2);
193   if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
194   *coordinates = ctx->coords;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "DMInterpolationGetVector"
200 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) {
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   PetscValidPointer(v, 2);
205   if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
206   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
207   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
208   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
209   ierr = VecSetFromOptions(*v);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "DMInterpolationRestoreVector"
215 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) {
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 PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
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) {
247       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) {
252         xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
253       }
254       for(comp = 0; comp < ctx->dof; ++comp) {
255         a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
256       }
257     }
258     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr);
259   }
260   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
261   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
262   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "QuadMap_Private"
268 PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
269 {
270   const PetscScalar*vertices = (const PetscScalar *) ctx;
271   const PetscScalar x0   = vertices[0];
272   const PetscScalar y0   = vertices[1];
273   const PetscScalar x1   = vertices[2];
274   const PetscScalar y1   = vertices[3];
275   const PetscScalar x2   = vertices[4];
276   const PetscScalar y2   = vertices[5];
277   const PetscScalar x3   = vertices[6];
278   const PetscScalar y3   = vertices[7];
279   const PetscScalar f_1  = x1 - x0;
280   const PetscScalar g_1  = y1 - y0;
281   const PetscScalar f_3  = x3 - x0;
282   const PetscScalar g_3  = y3 - y0;
283   const PetscScalar f_01 = x2 - x1 - x3 + x0;
284   const PetscScalar g_01 = y2 - y1 - y3 + y0;
285   PetscScalar      *ref, *real;
286   PetscErrorCode    ierr;
287 
288   PetscFunctionBegin;
289   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
290   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
291   {
292     const PetscScalar p0 = ref[0];
293     const PetscScalar p1 = ref[1];
294 
295     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
296     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
297   }
298   ierr = PetscLogFlops(28);CHKERRQ(ierr);
299   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
300   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "QuadJacobian_Private"
306 PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
307 {
308   const PetscScalar*vertices = (const PetscScalar *) ctx;
309   const PetscScalar x0   = vertices[0];
310   const PetscScalar y0   = vertices[1];
311   const PetscScalar x1   = vertices[2];
312   const PetscScalar y1   = vertices[3];
313   const PetscScalar x2   = vertices[4];
314   const PetscScalar y2   = vertices[5];
315   const PetscScalar x3   = vertices[6];
316   const PetscScalar y3   = vertices[7];
317   const PetscScalar f_01 = x2 - x1 - x3 + x0;
318   const PetscScalar g_01 = y2 - y1 - y3 + y0;
319   PetscScalar      *ref;
320   PetscErrorCode    ierr;
321 
322   PetscFunctionBegin;
323   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
324   {
325     const PetscScalar x = ref[0];
326     const PetscScalar y = ref[1];
327     const PetscInt    rows[2]   = {0, 1};
328     const PetscScalar values[4] = {(x1 - x0 + f_01*y) * 0.5, (x3 - x0 + f_01*x) * 0.5,
329                                    (y1 - y0 + g_01*y) * 0.5, (y3 - y0 + g_01*x) * 0.5};
330     ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
331   }
332   ierr = PetscLogFlops(30);CHKERRQ(ierr);
333   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
334   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
335   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
336   PetscFunctionReturn(0);
337 }
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "DMInterpolate_Quad_Private"
341 PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
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 = 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(dm, 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) {
395       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(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 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 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     const PetscScalar values[9] = {
534       (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0,
535       (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0,
536       (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0,
537       (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0,
538       (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0,
539       (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0,
540       (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0,
541       (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0,
542       (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0};
543     ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
544   }
545   ierr = PetscLogFlops(152);CHKERRQ(ierr);
546   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
547   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
548   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
549   PetscFunctionReturn(0);
550 }
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "DMInterpolate_Hex_Private"
554 PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
555   SNES           snes;
556   KSP            ksp;
557   PC             pc;
558   Vec            coordsLocal, r, ref, real;
559   Mat            J;
560   PetscScalar   *a, *coords;
561   PetscInt       p;
562   PetscErrorCode ierr;
563 
564   PetscFunctionBegin;
565   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
566   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
567   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
568   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
569   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
570   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
571   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
572   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
573   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
574   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
575   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
576   ierr = MatSetUp(J);CHKERRQ(ierr);
577   ierr = SNESSetFunction(snes, r, HexMap_Private, PETSC_NULL);CHKERRQ(ierr);
578   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, PETSC_NULL);CHKERRQ(ierr);
579   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
580   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
581   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
582   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
583 
584   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
585   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
586   for(p = 0; p < ctx->n; ++p) {
587     const PetscScalar *x, *vertices;
588     PetscScalar       *xi;
589     PetscReal          xir[3];
590     PetscInt           c = ctx->cells[p], comp, coordSize, xSize;
591 
592     /* Can make this do all points at once */
593     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
594     if (8*3 != coordSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);}
595     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
596     if (8*ctx->dof != xSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);}
597     ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
598     ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
599     ierr = VecGetArray(real, &xi);CHKERRQ(ierr);
600     xi[0] = coords[p*ctx->dim+0];
601     xi[1] = coords[p*ctx->dim+1];
602     xi[2] = coords[p*ctx->dim+2];
603     ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr);
604     ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr);
605     ierr = VecGetArray(ref, &xi);CHKERRQ(ierr);
606     xir[0] = PetscRealPart(xi[0]);
607     xir[1] = PetscRealPart(xi[1]);
608     xir[2] = PetscRealPart(xi[2]);
609     for (comp = 0; comp < ctx->dof; ++comp) {
610       a[p*ctx->dof+comp] =
611         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
612         x[1*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
613         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
614         x[3*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
615         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
616         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
617         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
618         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
619     }
620     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
621     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
622     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
623   }
624   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
625   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
626 
627   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
628   ierr = VecDestroy(&r);CHKERRQ(ierr);
629   ierr = VecDestroy(&ref);CHKERRQ(ierr);
630   ierr = VecDestroy(&real);CHKERRQ(ierr);
631   ierr = MatDestroy(&J);CHKERRQ(ierr);
632   PetscFunctionReturn(0);
633 }
634 
635 #undef __FUNCT__
636 #define __FUNCT__ "DMInterpolationEvaluate"
637 /*
638   Input Parameters:
639 + ctx - The DMInterpolationInfo context
640 . dm  - The DM
641 - x   - The local vector containing the field to be interpolated
642 
643   Output Parameters:
644 . v   - The vector containing the interpolated values
645 */
646 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) {
647   PetscInt       dim, coneSize, n;
648   PetscErrorCode ierr;
649 
650   PetscFunctionBegin;
651   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
652   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
653   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
654   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
655   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);}
656   if (n) {
657     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
658     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
659     if (dim == 2) {
660       if (coneSize == 3) {
661         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
662       } else if (coneSize == 4) {
663         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
664       } else {
665         SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
666       }
667     } else if (dim == 3) {
668       if (coneSize == 4) {
669         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
670       } else {
671         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
672       }
673     } else {
674       SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
675     }
676   }
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "DMInterpolationDestroy"
682 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) {
683   PetscErrorCode ierr;
684 
685   PetscFunctionBegin;
686   PetscValidPointer(ctx, 2);
687   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
688   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
689   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
690   ierr = PetscFree(*ctx);CHKERRQ(ierr);
691   *ctx = PETSC_NULL;
692   PetscFunctionReturn(0);
693 }
694