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