xref: /petsc/src/snes/utils/dmplexsnes.c (revision 26f7fa13d894e02c8192cf9a0d4955dba2e1288c)
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   PetscReal     *v0, *J, *invJ, detJ;
232   PetscScalar   *a, *coords;
233   PetscInt       p;
234   PetscErrorCode ierr;
235 
236   PetscFunctionBegin;
237   ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr);
238   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
239   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
240   for(p = 0; p < ctx->n; ++p) {
241     PetscInt           c = ctx->cells[p];
242     const PetscScalar *x;
243     PetscReal          xi[4];
244     PetscInt           d, f, comp;
245 
246     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
247     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
248     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr);
249     for(comp = 0; comp < ctx->dof; ++comp) {
250       a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
251     }
252     for(d = 0; d < ctx->dim; ++d) {
253       xi[d] = 0.0;
254       for(f = 0; f < ctx->dim; ++f) {
255         xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
256       }
257       for(comp = 0; comp < ctx->dof; ++comp) {
258         a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
259       }
260     }
261     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr);
262   }
263   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
264   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
265   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "QuadMap_Private"
271 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
272 {
273   const PetscScalar*vertices = (const PetscScalar *) ctx;
274   const PetscScalar x0   = vertices[0];
275   const PetscScalar y0   = vertices[1];
276   const PetscScalar x1   = vertices[2];
277   const PetscScalar y1   = vertices[3];
278   const PetscScalar x2   = vertices[4];
279   const PetscScalar y2   = vertices[5];
280   const PetscScalar x3   = vertices[6];
281   const PetscScalar y3   = vertices[7];
282   const PetscScalar f_1  = x1 - x0;
283   const PetscScalar g_1  = y1 - y0;
284   const PetscScalar f_3  = x3 - x0;
285   const PetscScalar g_3  = y3 - y0;
286   const PetscScalar f_01 = x2 - x1 - x3 + x0;
287   const PetscScalar g_01 = y2 - y1 - y3 + y0;
288   PetscScalar      *ref, *real;
289   PetscErrorCode    ierr;
290 
291   PetscFunctionBegin;
292   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
293   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
294   {
295     const PetscScalar p0 = ref[0];
296     const PetscScalar p1 = ref[1];
297 
298     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
299     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
300   }
301   ierr = PetscLogFlops(28);CHKERRQ(ierr);
302   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
303   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "QuadJacobian_Private"
309 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
310 {
311   const PetscScalar*vertices = (const PetscScalar *) ctx;
312   const PetscScalar x0   = vertices[0];
313   const PetscScalar y0   = vertices[1];
314   const PetscScalar x1   = vertices[2];
315   const PetscScalar y1   = vertices[3];
316   const PetscScalar x2   = vertices[4];
317   const PetscScalar y2   = vertices[5];
318   const PetscScalar x3   = vertices[6];
319   const PetscScalar y3   = vertices[7];
320   const PetscScalar f_01 = x2 - x1 - x3 + x0;
321   const PetscScalar g_01 = y2 - y1 - y3 + y0;
322   PetscScalar      *ref;
323   PetscErrorCode    ierr;
324 
325   PetscFunctionBegin;
326   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
327   {
328     const PetscScalar x = ref[0];
329     const PetscScalar y = ref[1];
330     const PetscInt    rows[2]   = {0, 1};
331     const PetscScalar values[4] = {(x1 - x0 + f_01*y) * 0.5, (x3 - x0 + f_01*x) * 0.5,
332                                    (y1 - y0 + g_01*y) * 0.5, (y3 - y0 + g_01*x) * 0.5};
333     ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
334   }
335   ierr = PetscLogFlops(30);CHKERRQ(ierr);
336   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
337   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
338   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339   PetscFunctionReturn(0);
340 }
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "DMInterpolate_Quad_Private"
344 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
345   DM             dmCoord;
346   SNES           snes;
347   KSP            ksp;
348   PC             pc;
349   Vec            coordsLocal, r, ref, real;
350   Mat            J;
351   PetscScalar   *a, *coords;
352   PetscInt       p;
353   PetscErrorCode ierr;
354 
355   PetscFunctionBegin;
356   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
357   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
358   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
359   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
360   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
361   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
362   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
363   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
364   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
365   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
366   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
367   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
368   ierr = MatSetUp(J);CHKERRQ(ierr);
369   ierr = SNESSetFunction(snes, r, QuadMap_Private, PETSC_NULL);CHKERRQ(ierr);
370   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, PETSC_NULL);CHKERRQ(ierr);
371   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
372   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
373   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
374   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
375 
376   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
377   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
378   for (p = 0; p < ctx->n; ++p) {
379     const PetscScalar *x, *vertices;
380     PetscScalar       *xi;
381     PetscReal          xir[2];
382     PetscInt           c = ctx->cells[p], comp, coordSize, xSize;
383 
384     /* Can make this do all points at once */
385     ierr = DMPlexVecGetClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
386     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
387     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
388     if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);
389     ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
390     ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
391     ierr = VecGetArray(real, &xi);CHKERRQ(ierr);
392     xi[0] = coords[p*ctx->dim+0];
393     xi[1] = coords[p*ctx->dim+1];
394     ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr);
395     ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr);
396     ierr = VecGetArray(ref, &xi);CHKERRQ(ierr);
397     xir[0] = PetscRealPart(xi[0]);
398     xir[1] = PetscRealPart(xi[1]);
399     for (comp = 0; comp < ctx->dof; ++comp) {
400       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];
401     }
402     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
403     ierr = DMPlexVecRestoreClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
404     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
405   }
406   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
407   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
408 
409   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
410   ierr = VecDestroy(&r);CHKERRQ(ierr);
411   ierr = VecDestroy(&ref);CHKERRQ(ierr);
412   ierr = VecDestroy(&real);CHKERRQ(ierr);
413   ierr = MatDestroy(&J);CHKERRQ(ierr);
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "HexMap_Private"
419 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
420 {
421   const PetscScalar*vertices = (const PetscScalar *) ctx;
422   const PetscScalar x0 = vertices[0];
423   const PetscScalar y0 = vertices[1];
424   const PetscScalar z0 = vertices[2];
425   const PetscScalar x1 = vertices[3];
426   const PetscScalar y1 = vertices[4];
427   const PetscScalar z1 = vertices[5];
428   const PetscScalar x2 = vertices[6];
429   const PetscScalar y2 = vertices[7];
430   const PetscScalar z2 = vertices[8];
431   const PetscScalar x3 = vertices[9];
432   const PetscScalar y3 = vertices[10];
433   const PetscScalar z3 = vertices[11];
434   const PetscScalar x4 = vertices[12];
435   const PetscScalar y4 = vertices[13];
436   const PetscScalar z4 = vertices[14];
437   const PetscScalar x5 = vertices[15];
438   const PetscScalar y5 = vertices[16];
439   const PetscScalar z5 = vertices[17];
440   const PetscScalar x6 = vertices[18];
441   const PetscScalar y6 = vertices[19];
442   const PetscScalar z6 = vertices[20];
443   const PetscScalar x7 = vertices[21];
444   const PetscScalar y7 = vertices[22];
445   const PetscScalar z7 = vertices[23];
446   const PetscScalar f_1 = x1 - x0;
447   const PetscScalar g_1 = y1 - y0;
448   const PetscScalar h_1 = z1 - z0;
449   const PetscScalar f_3 = x3 - x0;
450   const PetscScalar g_3 = y3 - y0;
451   const PetscScalar h_3 = z3 - z0;
452   const PetscScalar f_4 = x4 - x0;
453   const PetscScalar g_4 = y4 - y0;
454   const PetscScalar h_4 = z4 - z0;
455   const PetscScalar f_01 = x2 - x1 - x3 + x0;
456   const PetscScalar g_01 = y2 - y1 - y3 + y0;
457   const PetscScalar h_01 = z2 - z1 - z3 + z0;
458   const PetscScalar f_12 = x7 - x3 - x4 + x0;
459   const PetscScalar g_12 = y7 - y3 - y4 + y0;
460   const PetscScalar h_12 = z7 - z3 - z4 + z0;
461   const PetscScalar f_02 = x5 - x1 - x4 + x0;
462   const PetscScalar g_02 = y5 - y1 - y4 + y0;
463   const PetscScalar h_02 = z5 - z1 - z4 + z0;
464   const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
465   const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
466   const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
467   PetscScalar      *ref, *real;
468   PetscErrorCode    ierr;
469 
470   PetscFunctionBegin;
471   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
472   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
473   {
474     const PetscScalar p0 = ref[0];
475     const PetscScalar p1 = ref[1];
476     const PetscScalar p2 = ref[2];
477 
478     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;
479     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;
480     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;
481   }
482   ierr = PetscLogFlops(114);CHKERRQ(ierr);
483   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
484   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "HexJacobian_Private"
490 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
491 {
492   const PetscScalar*vertices = (const PetscScalar *) ctx;
493   const PetscScalar x0 = vertices[0];
494   const PetscScalar y0 = vertices[1];
495   const PetscScalar z0 = vertices[2];
496   const PetscScalar x1 = vertices[3];
497   const PetscScalar y1 = vertices[4];
498   const PetscScalar z1 = vertices[5];
499   const PetscScalar x2 = vertices[6];
500   const PetscScalar y2 = vertices[7];
501   const PetscScalar z2 = vertices[8];
502   const PetscScalar x3 = vertices[9];
503   const PetscScalar y3 = vertices[10];
504   const PetscScalar z3 = vertices[11];
505   const PetscScalar x4 = vertices[12];
506   const PetscScalar y4 = vertices[13];
507   const PetscScalar z4 = vertices[14];
508   const PetscScalar x5 = vertices[15];
509   const PetscScalar y5 = vertices[16];
510   const PetscScalar z5 = vertices[17];
511   const PetscScalar x6 = vertices[18];
512   const PetscScalar y6 = vertices[19];
513   const PetscScalar z6 = vertices[20];
514   const PetscScalar x7 = vertices[21];
515   const PetscScalar y7 = vertices[22];
516   const PetscScalar z7 = vertices[23];
517   const PetscScalar f_xy = x2 - x1 - x3 + x0;
518   const PetscScalar g_xy = y2 - y1 - y3 + y0;
519   const PetscScalar h_xy = z2 - z1 - z3 + z0;
520   const PetscScalar f_yz = x7 - x3 - x4 + x0;
521   const PetscScalar g_yz = y7 - y3 - y4 + y0;
522   const PetscScalar h_yz = z7 - z3 - z4 + z0;
523   const PetscScalar f_xz = x5 - x1 - x4 + x0;
524   const PetscScalar g_xz = y5 - y1 - y4 + y0;
525   const PetscScalar h_xz = z5 - z1 - z4 + z0;
526   const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
527   const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
528   const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
529   PetscScalar      *ref;
530   PetscErrorCode    ierr;
531 
532   PetscFunctionBegin;
533   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
534   {
535     const PetscScalar x = ref[0];
536     const PetscScalar y = ref[1];
537     const PetscScalar z = ref[2];
538     const PetscInt    rows[3]   = {0, 1, 2};
539     const PetscScalar values[9] = {
540       (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0,
541       (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0,
542       (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0,
543       (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0,
544       (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0,
545       (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0,
546       (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0,
547       (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0,
548       (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0};
549     ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
550   }
551   ierr = PetscLogFlops(152);CHKERRQ(ierr);
552   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
553   ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
554   ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
555   PetscFunctionReturn(0);
556 }
557 
558 #undef __FUNCT__
559 #define __FUNCT__ "DMInterpolate_Hex_Private"
560 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) {
561   DM             dmCoord;
562   SNES           snes;
563   KSP            ksp;
564   PC             pc;
565   Vec            coordsLocal, r, ref, real;
566   Mat            J;
567   PetscScalar   *a, *coords;
568   PetscInt       p;
569   PetscErrorCode ierr;
570 
571   PetscFunctionBegin;
572   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
573   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
574   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
575   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
576   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
577   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
578   ierr = VecSetFromOptions(r);CHKERRQ(ierr);
579   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
580   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
581   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
582   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
583   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
584   ierr = MatSetUp(J);CHKERRQ(ierr);
585   ierr = SNESSetFunction(snes, r, HexMap_Private, PETSC_NULL);CHKERRQ(ierr);
586   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, PETSC_NULL);CHKERRQ(ierr);
587   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
588   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
589   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
590   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
591 
592   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
593   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
594   for(p = 0; p < ctx->n; ++p) {
595     const PetscScalar *x, *vertices;
596     PetscScalar       *xi;
597     PetscReal          xir[3];
598     PetscInt           c = ctx->cells[p], comp, coordSize, xSize;
599 
600     /* Can make this do all points at once */
601     ierr = DMPlexVecGetClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
602     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
603     ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
604     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
605     ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
606     ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void *) vertices);CHKERRQ(ierr);
607     ierr = VecGetArray(real, &xi);CHKERRQ(ierr);
608     xi[0] = coords[p*ctx->dim+0];
609     xi[1] = coords[p*ctx->dim+1];
610     xi[2] = coords[p*ctx->dim+2];
611     ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr);
612     ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr);
613     ierr = VecGetArray(ref, &xi);CHKERRQ(ierr);
614     xir[0] = PetscRealPart(xi[0]);
615     xir[1] = PetscRealPart(xi[1]);
616     xir[2] = PetscRealPart(xi[2]);
617     for (comp = 0; comp < ctx->dof; ++comp) {
618       a[p*ctx->dof+comp] =
619         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
620         x[1*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
621         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
622         x[3*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
623         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
624         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
625         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
626         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
627     }
628     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
629     ierr = DMPlexVecRestoreClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
630     ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
631   }
632   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
633   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
634 
635   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
636   ierr = VecDestroy(&r);CHKERRQ(ierr);
637   ierr = VecDestroy(&ref);CHKERRQ(ierr);
638   ierr = VecDestroy(&real);CHKERRQ(ierr);
639   ierr = MatDestroy(&J);CHKERRQ(ierr);
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "DMInterpolationEvaluate"
645 /*
646   Input Parameters:
647 + ctx - The DMInterpolationInfo context
648 . dm  - The DM
649 - x   - The local vector containing the field to be interpolated
650 
651   Output Parameters:
652 . v   - The vector containing the interpolated values
653 */
654 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
655 {
656   PetscInt       dim, coneSize, n;
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
661   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
662   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
663   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
664   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);
665   if (n) {
666     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
667     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
668     if (dim == 2) {
669       if (coneSize == 3) {
670         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
671       } else if (coneSize == 4) {
672         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
673       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
674     } else if (dim == 3) {
675       if (coneSize == 4) {
676         ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr);
677       } else {
678         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
679       }
680     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
681   }
682   PetscFunctionReturn(0);
683 }
684 
685 #undef __FUNCT__
686 #define __FUNCT__ "DMInterpolationDestroy"
687 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
688 {
689   PetscErrorCode ierr;
690 
691   PetscFunctionBegin;
692   PetscValidPointer(ctx, 2);
693   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
694   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
695   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
696   ierr = PetscFree(*ctx);CHKERRQ(ierr);
697   *ctx = PETSC_NULL;
698   PetscFunctionReturn(0);
699 }
700