Actual source code: plexgeometry.c
petsc-3.8.4 2018-03-24
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/petscfeimpl.h>
3: #include <petscblaslapack.h>
4: #include <petsctime.h>
6: static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection)
7: {
8: const PetscReal p0_x = segmentA[0*2+0];
9: const PetscReal p0_y = segmentA[0*2+1];
10: const PetscReal p1_x = segmentA[1*2+0];
11: const PetscReal p1_y = segmentA[1*2+1];
12: const PetscReal p2_x = segmentB[0*2+0];
13: const PetscReal p2_y = segmentB[0*2+1];
14: const PetscReal p3_x = segmentB[1*2+0];
15: const PetscReal p3_y = segmentB[1*2+1];
16: const PetscReal s1_x = p1_x - p0_x;
17: const PetscReal s1_y = p1_y - p0_y;
18: const PetscReal s2_x = p3_x - p2_x;
19: const PetscReal s2_y = p3_y - p2_y;
20: const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y);
23: *hasIntersection = PETSC_FALSE;
24: /* Non-parallel lines */
25: if (denom != 0.0) {
26: const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom;
27: const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom;
29: if (s >= 0 && s <= 1 && t >= 0 && t <= 1) {
30: *hasIntersection = PETSC_TRUE;
31: if (intersection) {
32: intersection[0] = p0_x + (t * s1_x);
33: intersection[1] = p0_y + (t * s1_y);
34: }
35: }
36: }
37: return(0);
38: }
40: static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
41: {
42: const PetscInt embedDim = 2;
43: const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON;
44: PetscReal x = PetscRealPart(point[0]);
45: PetscReal y = PetscRealPart(point[1]);
46: PetscReal v0[2], J[4], invJ[4], detJ;
47: PetscReal xi, eta;
48: PetscErrorCode ierr;
51: DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
52: xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
53: eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
55: if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c;
56: else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
57: return(0);
58: }
60: static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[])
61: {
62: const PetscInt embedDim = 2;
63: PetscReal x = PetscRealPart(point[0]);
64: PetscReal y = PetscRealPart(point[1]);
65: PetscReal v0[2], J[4], invJ[4], detJ;
66: PetscReal xi, eta, r;
67: PetscErrorCode ierr;
70: DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
71: xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]);
72: eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]);
74: xi = PetscMax(xi, 0.0);
75: eta = PetscMax(eta, 0.0);
76: if (xi + eta > 2.0) {
77: r = (xi + eta)/2.0;
78: xi /= r;
79: eta /= r;
80: }
81: cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0];
82: cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1];
83: return(0);
84: }
86: static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
87: {
88: PetscSection coordSection;
89: Vec coordsLocal;
90: PetscScalar *coords = NULL;
91: const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0};
92: PetscReal x = PetscRealPart(point[0]);
93: PetscReal y = PetscRealPart(point[1]);
94: PetscInt crossings = 0, f;
95: PetscErrorCode ierr;
98: DMGetCoordinatesLocal(dm, &coordsLocal);
99: DMGetCoordinateSection(dm, &coordSection);
100: DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);
101: for (f = 0; f < 4; ++f) {
102: PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]);
103: PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]);
104: PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]);
105: PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]);
106: PetscReal slope = (y_j - y_i) / (x_j - x_i);
107: PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE;
108: PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE;
109: PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE;
110: if ((cond1 || cond2) && above) ++crossings;
111: }
112: if (crossings % 2) *cell = c;
113: else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
114: DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);
115: return(0);
116: }
118: static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
119: {
120: const PetscInt embedDim = 3;
121: PetscReal v0[3], J[9], invJ[9], detJ;
122: PetscReal x = PetscRealPart(point[0]);
123: PetscReal y = PetscRealPart(point[1]);
124: PetscReal z = PetscRealPart(point[2]);
125: PetscReal xi, eta, zeta;
129: DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
130: xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]);
131: eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]);
132: zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]);
134: if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c;
135: else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
136: return(0);
137: }
139: static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell)
140: {
141: PetscSection coordSection;
142: Vec coordsLocal;
143: PetscScalar *coords;
144: const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5,
145: 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4};
146: PetscBool found = PETSC_TRUE;
147: PetscInt f;
151: DMGetCoordinatesLocal(dm, &coordsLocal);
152: DMGetCoordinateSection(dm, &coordSection);
153: DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);
154: for (f = 0; f < 6; ++f) {
155: /* Check the point is under plane */
156: /* Get face normal */
157: PetscReal v_i[3];
158: PetscReal v_j[3];
159: PetscReal normal[3];
160: PetscReal pp[3];
161: PetscReal dot;
163: v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]);
164: v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]);
165: v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]);
166: v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]);
167: v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]);
168: v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]);
169: normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1];
170: normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2];
171: normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0];
172: pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]);
173: pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]);
174: pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]);
175: dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2];
177: /* Check that projected point is in face (2D location problem) */
178: if (dot < 0.0) {
179: found = PETSC_FALSE;
180: break;
181: }
182: }
183: if (found) *cell = c;
184: else *cell = DMLOCATEPOINT_POINT_NOT_FOUND;
185: DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);
186: return(0);
187: }
189: static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[])
190: {
191: PetscInt d;
194: box->dim = dim;
195: for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]);
196: return(0);
197: }
199: PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box)
200: {
204: PetscMalloc1(1, box);
205: PetscGridHashInitialize_Internal(*box, dim, point);
206: return(0);
207: }
209: PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[])
210: {
211: PetscInt d;
214: for (d = 0; d < box->dim; ++d) {
215: box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d]));
216: box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d]));
217: }
218: return(0);
219: }
221: /*
222: PetscGridHashSetGrid - Divide the grid into boxes
224: Not collective
226: Input Parameters:
227: + box - The grid hash object
228: . n - The number of boxes in each dimension, or PETSC_DETERMINE
229: - h - The box size in each dimension, only used if n[d] == PETSC_DETERMINE
231: Level: developer
233: .seealso: PetscGridHashCreate()
234: */
235: PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[])
236: {
237: PetscInt d;
240: for (d = 0; d < box->dim; ++d) {
241: box->extent[d] = box->upper[d] - box->lower[d];
242: if (n[d] == PETSC_DETERMINE) {
243: box->h[d] = h[d];
244: box->n[d] = PetscCeilReal(box->extent[d]/h[d]);
245: } else {
246: box->n[d] = n[d];
247: box->h[d] = box->extent[d]/n[d];
248: }
249: }
250: return(0);
251: }
253: /*
254: PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point
256: Not collective
258: Input Parameters:
259: + box - The grid hash object
260: . numPoints - The number of input points
261: - points - The input point coordinates
263: Output Parameters:
264: + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
265: - boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL
267: Level: developer
269: .seealso: PetscGridHashCreate()
270: */
271: PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[])
272: {
273: const PetscReal *lower = box->lower;
274: const PetscReal *upper = box->upper;
275: const PetscReal *h = box->h;
276: const PetscInt *n = box->n;
277: const PetscInt dim = box->dim;
278: PetscInt d, p;
281: for (p = 0; p < numPoints; ++p) {
282: for (d = 0; d < dim; ++d) {
283: PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
285: if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
286: if (dbox == -1 && PetscAbsReal(PetscRealPart(points[p*dim+d]) - lower[d]) < 1.0e-9) dbox = 0;
287: if (dbox < 0 || dbox >= n[d]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %d (%g, %g, %g) is outside of our bounding box",
288: p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0);
289: dboxes[p*dim+d] = dbox;
290: }
291: if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
292: }
293: return(0);
294: }
296: /*
297: PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point
298:
299: Not collective
300:
301: Input Parameters:
302: + box - The grid hash object
303: . numPoints - The number of input points
304: - points - The input point coordinates
305:
306: Output Parameters:
307: + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim)
308: . boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL
309: - found - Flag indicating if point was located within a box
310:
311: Level: developer
312:
313: .seealso: PetscGridHashGetEnclosingBox()
314: */
315: PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found)
316: {
317: const PetscReal *lower = box->lower;
318: const PetscReal *upper = box->upper;
319: const PetscReal *h = box->h;
320: const PetscInt *n = box->n;
321: const PetscInt dim = box->dim;
322: PetscInt d, p;
323:
325: *found = PETSC_FALSE;
326: for (p = 0; p < numPoints; ++p) {
327: for (d = 0; d < dim; ++d) {
328: PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]);
329:
330: if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1;
331: if (dbox < 0 || dbox >= n[d]) {
332: return(0);
333: }
334: dboxes[p*dim+d] = dbox;
335: }
336: if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1];
337: }
338: *found = PETSC_TRUE;
339: return(0);
340: }
342: PetscErrorCode PetscGridHashDestroy(PetscGridHash *box)
343: {
347: if (*box) {
348: PetscSectionDestroy(&(*box)->cellSection);
349: ISDestroy(&(*box)->cells);
350: DMLabelDestroy(&(*box)->cellsSparse);
351: }
352: PetscFree(*box);
353: return(0);
354: }
356: PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell)
357: {
358: PetscInt coneSize;
362: switch (dim) {
363: case 2:
364: DMPlexGetConeSize(dm, cellStart, &coneSize);
365: switch (coneSize) {
366: case 3:
367: DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);
368: break;
369: case 4:
370: DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);
371: break;
372: default:
373: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
374: }
375: break;
376: case 3:
377: DMPlexGetConeSize(dm, cellStart, &coneSize);
378: switch (coneSize) {
379: case 4:
380: DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);
381: break;
382: case 6:
383: DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);
384: break;
385: default:
386: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize);
387: }
388: break;
389: default:
390: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim);
391: }
392: return(0);
393: }
395: /*
396: DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point
397: */
398: PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[])
399: {
400: PetscInt coneSize;
404: switch (dim) {
405: case 2:
406: DMPlexGetConeSize(dm, cell, &coneSize);
407: switch (coneSize) {
408: case 3:
409: DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);
410: break;
411: #if 0
412: case 4:
413: DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);
414: break;
415: #endif
416: default:
417: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
418: }
419: break;
420: #if 0
421: case 3:
422: DMPlexGetConeSize(dm, cell, &coneSize);
423: switch (coneSize) {
424: case 4:
425: DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);
426: break;
427: case 6:
428: DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);
429: break;
430: default:
431: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize);
432: }
433: break;
434: #endif
435: default:
436: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim);
437: }
438: return(0);
439: }
441: /*
442: DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex
444: Collective on DM
446: Input Parameter:
447: . dm - The Plex
449: Output Parameter:
450: . localBox - The grid hash object
452: Level: developer
454: .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox()
455: */
456: PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox)
457: {
458: MPI_Comm comm;
459: PetscGridHash lbox;
460: Vec coordinates;
461: PetscSection coordSection;
462: Vec coordsLocal;
463: const PetscScalar *coords;
464: PetscInt *dboxes, *boxes;
465: PetscInt n[3] = {10, 10, 10};
466: PetscInt dim, N, cStart, cEnd, cMax, c, i;
467: PetscErrorCode ierr;
470: PetscObjectGetComm((PetscObject) dm, &comm);
471: DMGetCoordinatesLocal(dm, &coordinates);
472: DMGetCoordinateDim(dm, &dim);
473: if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D");
474: VecGetLocalSize(coordinates, &N);
475: VecGetArrayRead(coordinates, &coords);
476: PetscGridHashCreate(comm, dim, coords, &lbox);
477: for (i = 0; i < N; i += dim) {PetscGridHashEnlarge(lbox, &coords[i]);}
478: VecRestoreArrayRead(coordinates, &coords);
479: PetscOptionsGetInt(NULL,NULL,"-dm_plex_hash_box_nijk",&n[0],NULL);
480: n[1] = n[0];
481: n[2] = n[0];
482: PetscGridHashSetGrid(lbox, n, NULL);
483: #if 0
484: /* Could define a custom reduction to merge these */
485: MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);
486: MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);
487: #endif
488: /* Is there a reason to snap the local bounding box to a division of the global box? */
489: /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */
490: /* Create label */
491: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
492: DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
493: if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
494: DMLabelCreate("cells", &lbox->cellsSparse);
495: DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);
496: /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */
497: DMGetCoordinatesLocal(dm, &coordsLocal);
498: DMGetCoordinateSection(dm, &coordSection);
499: PetscCalloc2(16 * dim, &dboxes, 16, &boxes);
500: for (c = cStart; c < cEnd; ++c) {
501: const PetscReal *h = lbox->h;
502: PetscScalar *ccoords = NULL;
503: PetscInt csize = 0;
504: PetscScalar point[3];
505: PetscInt dlim[6], d, e, i, j, k;
507: /* Find boxes enclosing each vertex */
508: DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);
509: PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);
510: /* Mark cells containing the vertices */
511: for (e = 0; e < csize/dim; ++e) {DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);}
512: /* Get grid of boxes containing these */
513: for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];}
514: for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;}
515: for (e = 1; e < dim+1; ++e) {
516: for (d = 0; d < dim; ++d) {
517: dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]);
518: dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]);
519: }
520: }
521: /* Check for intersection of box with cell */
522: for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) {
523: for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) {
524: for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) {
525: const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i;
526: PetscScalar cpoint[3];
527: PetscInt cell, edge, ii, jj, kk;
529: /* Check whether cell contains any vertex of these subboxes TODO vectorize this */
530: for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) {
531: for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) {
532: for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) {
534: DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);
535: if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box); ii = jj = kk = 2;}
536: }
537: }
538: }
539: /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */
540: for (edge = 0; edge < dim+1; ++edge) {
541: PetscReal segA[6], segB[6];
543: if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected dim %d > 3",dim);
544: for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);}
545: for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) {
546: if (dim > 2) {segB[2] = PetscRealPart(point[2]);
547: segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];}
548: for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) {
549: if (dim > 1) {segB[1] = PetscRealPart(point[1]);
550: segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];}
551: for (ii = 0; ii < 2; ++ii) {
552: PetscBool intersects;
554: segB[0] = PetscRealPart(point[0]);
555: segB[dim+0] = PetscRealPart(point[0]) + ii*h[0];
556: DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);
557: if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box); edge = ii = jj = kk = dim+1;}
558: }
559: }
560: }
561: }
562: }
563: }
564: }
565: DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);
566: }
567: PetscFree2(dboxes, boxes);
568: DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);
569: DMLabelDestroy(&lbox->cellsSparse);
570: *localBox = lbox;
571: return(0);
572: }
574: PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF)
575: {
576: DM_Plex *mesh = (DM_Plex *) dm->data;
577: PetscBool hash = mesh->useHashLocation, reuse = PETSC_FALSE;
578: PetscInt bs, numPoints, p, numFound, *found = NULL;
579: PetscInt dim, cStart, cEnd, cMax, numCells, c, d;
580: const PetscInt *boxCells;
581: PetscSFNode *cells;
582: PetscScalar *a;
583: PetscMPIInt result;
584: PetscLogDouble t0,t1;
585: PetscReal gmin[3],gmax[3];
586: PetscInt terminating_query_type[] = { 0, 0, 0 };
587: PetscErrorCode ierr;
590: PetscTime(&t0);
591: if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it.");
592: DMGetCoordinateDim(dm, &dim);
593: VecGetBlockSize(v, &bs);
594: MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);
595: if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported");
596: if (bs != dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Block size for point vector %D must be the mesh coordinate dimension %D", bs, dim);
597: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
598: DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
599: if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
600: VecGetLocalSize(v, &numPoints);
601: VecGetArray(v, &a);
602: numPoints /= bs;
603: {
604: const PetscSFNode *sf_cells;
605:
606: PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);
607: if (sf_cells) {
608: PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");
609: cells = (PetscSFNode*)sf_cells;
610: reuse = PETSC_TRUE;
611: } else {
612: PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");
613: PetscMalloc1(numPoints, &cells);
614: /* initialize cells if created */
615: for (p=0; p<numPoints; p++) {
616: cells[p].rank = 0;
617: cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
618: }
619: }
620: }
621: /* define domain bounding box */
622: {
623: Vec coorglobal;
624:
625: DMGetCoordinates(dm,&coorglobal);
626: VecStrideMaxAll(coorglobal,NULL,gmax);
627: VecStrideMinAll(coorglobal,NULL,gmin);
628: }
629: if (hash) {
630: if (!mesh->lbox) {PetscInfo(dm, "Initializing grid hashing");DMPlexComputeGridHash_Internal(dm, &mesh->lbox);}
631: /* Designate the local box for each point */
632: /* Send points to correct process */
633: /* Search cells that lie in each subbox */
634: /* Should we bin points before doing search? */
635: ISGetIndices(mesh->lbox->cells, &boxCells);
636: }
637: for (p = 0, numFound = 0; p < numPoints; ++p) {
638: const PetscScalar *point = &a[p*bs];
639: PetscInt dbin[3] = {-1,-1,-1}, bin, cell = -1, cellOffset;
640: PetscBool point_outside_domain = PETSC_FALSE;
642: /* check bounding box of domain */
643: for (d=0; d<dim; d++) {
644: if (PetscRealPart(point[d]) < gmin[d]) { point_outside_domain = PETSC_TRUE; break; }
645: if (PetscRealPart(point[d]) > gmax[d]) { point_outside_domain = PETSC_TRUE; break; }
646: }
647: if (point_outside_domain) {
648: cells[p].rank = 0;
649: cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
650: terminating_query_type[0]++;
651: continue;
652: }
653:
654: /* check initial values in cells[].index - abort early if found */
655: if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
656: c = cells[p].index;
657: cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
658: DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);
659: if (cell >= 0) {
660: cells[p].rank = 0;
661: cells[p].index = cell;
662: numFound++;
663: }
664: }
665: if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
666: terminating_query_type[1]++;
667: continue;
668: }
669:
670: if (hash) {
671: PetscBool found_box;
672:
673: /* allow for case that point is outside box - abort early */
674: PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);
675: if (found_box) {
676: /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */
677: PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);
678: PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);
679: for (c = cellOffset; c < cellOffset + numCells; ++c) {
680: DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);
681: if (cell >= 0) {
682: cells[p].rank = 0;
683: cells[p].index = cell;
684: numFound++;
685: terminating_query_type[2]++;
686: break;
687: }
688: }
689: }
690: } else {
691: for (c = cStart; c < cEnd; ++c) {
692: DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);
693: if (cell >= 0) {
694: cells[p].rank = 0;
695: cells[p].index = cell;
696: numFound++;
697: terminating_query_type[2]++;
698: break;
699: }
700: }
701: }
702: }
703: if (hash) {ISRestoreIndices(mesh->lbox->cells, &boxCells);}
704: if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) {
705: for (p = 0; p < numPoints; p++) {
706: const PetscScalar *point = &a[p*bs];
707: PetscReal cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL;
708: PetscInt dbin[3] = {-1,-1,-1}, bin, cellOffset, d;
710: if (cells[p].index < 0) {
711: ++numFound;
712: PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);
713: PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);
714: PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);
715: for (c = cellOffset; c < cellOffset + numCells; ++c) {
716: DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);
717: for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]);
718: dist = DMPlex_NormD_Internal(dim, diff);
719: if (dist < distMax) {
720: for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d];
721: cells[p].rank = 0;
722: cells[p].index = boxCells[c];
723: distMax = dist;
724: break;
725: }
726: }
727: }
728: }
729: }
730: /* This code is only be relevant when interfaced to parallel point location */
731: /* Check for highest numbered proc that claims a point (do we care?) */
732: if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) {
733: PetscMalloc1(numFound,&found);
734: for (p = 0, numFound = 0; p < numPoints; p++) {
735: if (cells[p].rank >= 0 && cells[p].index >= 0) {
736: if (numFound < p) {
737: cells[numFound] = cells[p];
738: }
739: found[numFound++] = p;
740: }
741: }
742: }
743: VecRestoreArray(v, &a);
744: if (!reuse) {
745: PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);
746: }
747: PetscTime(&t1);
748: if (hash) {
749: PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [hash]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);
750: } else {
751: PetscInfo3(dm,"[DMLocatePoints_Plex] terminating_query_type : %D [outside domain] : %D [inside intial cell] : %D [brute-force]\n",terminating_query_type[0],terminating_query_type[1],terminating_query_type[2]);
752: }
753: PetscInfo3(dm,"[DMLocatePoints_Plex] npoints %D : time(rank0) %1.2e (sec): points/sec %1.4e\n",numPoints,t1-t0,(double)((double)numPoints/(t1-t0)));
754: return(0);
755: }
757: /*@C
758: DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates
760: Not collective
762: Input Parameter:
763: . coords - The coordinates of a segment
765: Output Parameters:
766: + coords - The new y-coordinate, and 0 for x
767: - R - The rotation which accomplishes the projection
769: Level: developer
771: .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D()
772: @*/
773: PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[])
774: {
775: const PetscReal x = PetscRealPart(coords[2] - coords[0]);
776: const PetscReal y = PetscRealPart(coords[3] - coords[1]);
777: const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r;
780: R[0] = c; R[1] = -s;
781: R[2] = s; R[3] = c;
782: coords[0] = 0.0;
783: coords[1] = r;
784: return(0);
785: }
787: /*@C
788: DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates
790: Not collective
792: Input Parameter:
793: . coords - The coordinates of a segment
795: Output Parameters:
796: + coords - The new y-coordinate, and 0 for x and z
797: - R - The rotation which accomplishes the projection
799: Note: This uses the basis completion described by Frisvad in http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html, DOI:10.1080/2165347X.2012.689606
801: Level: developer
803: .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D()
804: @*/
805: PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
806: {
807: PetscReal x = PetscRealPart(coords[3] - coords[0]);
808: PetscReal y = PetscRealPart(coords[4] - coords[1]);
809: PetscReal z = PetscRealPart(coords[5] - coords[2]);
810: PetscReal r = PetscSqrtReal(x*x + y*y + z*z);
811: PetscReal rinv = 1. / r;
814: x *= rinv; y *= rinv; z *= rinv;
815: if (x > 0.) {
816: PetscReal inv1pX = 1./ (1. + x);
818: R[0] = x; R[1] = -y; R[2] = -z;
819: R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX;
820: R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX;
821: }
822: else {
823: PetscReal inv1mX = 1./ (1. - x);
825: R[0] = x; R[1] = z; R[2] = y;
826: R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX;
827: R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX;
828: }
829: coords[0] = 0.0;
830: coords[1] = r;
831: return(0);
832: }
834: /*@
835: DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates
837: Not collective
839: Input Parameter:
840: . coords - The coordinates of a segment
842: Output Parameters:
843: + coords - The new y- and z-coordinates, and 0 for x
844: - R - The rotation which accomplishes the projection
846: Level: developer
848: .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()
849: @*/
850: PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])
851: {
852: PetscReal x1[3], x2[3], n[3], norm;
853: PetscReal x1p[3], x2p[3], xnp[3];
854: PetscReal sqrtz, alpha;
855: const PetscInt dim = 3;
856: PetscInt d, e, p;
859: /* 0) Calculate normal vector */
860: for (d = 0; d < dim; ++d) {
861: x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]);
862: x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]);
863: }
864: n[0] = x1[1]*x2[2] - x1[2]*x2[1];
865: n[1] = x1[2]*x2[0] - x1[0]*x2[2];
866: n[2] = x1[0]*x2[1] - x1[1]*x2[0];
867: norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]);
868: n[0] /= norm;
869: n[1] /= norm;
870: n[2] /= norm;
871: /* 1) Take the normal vector and rotate until it is \hat z
873: Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then
875: R = / alpha nx nz alpha ny nz -1/alpha \
876: | -alpha ny alpha nx 0 |
877: \ nx ny nz /
879: will rotate the normal vector to \hat z
880: */
881: sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]);
882: /* Check for n = z */
883: if (sqrtz < 1.0e-10) {
884: const PetscInt s = PetscSign(n[2]);
885: /* If nz < 0, rotate 180 degrees around x-axis */
886: for (p = 3; p < coordSize/3; ++p) {
887: coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]);
888: coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s;
889: }
890: coords[0] = 0.0;
891: coords[1] = 0.0;
892: coords[2] = x1[0];
893: coords[3] = x1[1] * s;
894: coords[4] = x2[0];
895: coords[5] = x2[1] * s;
896: R[0] = 1.0; R[1] = 0.0; R[2] = 0.0;
897: R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0;
898: R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s;
899: return(0);
900: }
901: alpha = 1.0/sqrtz;
902: R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz;
903: R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0;
904: R[6] = n[0]; R[7] = n[1]; R[8] = n[2];
905: for (d = 0; d < dim; ++d) {
906: x1p[d] = 0.0;
907: x2p[d] = 0.0;
908: for (e = 0; e < dim; ++e) {
909: x1p[d] += R[d*dim+e]*x1[e];
910: x2p[d] += R[d*dim+e]*x2[e];
911: }
912: }
913: if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
914: if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated");
915: /* 2) Project to (x, y) */
916: for (p = 3; p < coordSize/3; ++p) {
917: for (d = 0; d < dim; ++d) {
918: xnp[d] = 0.0;
919: for (e = 0; e < dim; ++e) {
920: xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]);
921: }
922: if (d < dim-1) coords[p*2+d] = xnp[d];
923: }
924: }
925: coords[0] = 0.0;
926: coords[1] = 0.0;
927: coords[2] = x1p[0];
928: coords[3] = x1p[1];
929: coords[4] = x2p[0];
930: coords[5] = x2p[1];
931: /* Output R^T which rotates \hat z to the input normal */
932: for (d = 0; d < dim; ++d) {
933: for (e = d+1; e < dim; ++e) {
934: PetscReal tmp;
936: tmp = R[d*dim+e];
937: R[d*dim+e] = R[e*dim+d];
938: R[e*dim+d] = tmp;
939: }
940: }
941: return(0);
942: }
944: PETSC_UNUSED
945: PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[])
946: {
947: /* Signed volume is 1/2 the determinant
949: | 1 1 1 |
950: | x0 x1 x2 |
951: | y0 y1 y2 |
953: but if x0,y0 is the origin, we have
955: | x1 x2 |
956: | y1 y2 |
957: */
958: const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1];
959: const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1];
960: PetscReal M[4], detM;
961: M[0] = x1; M[1] = x2;
962: M[2] = y1; M[3] = y2;
963: DMPlex_Det2D_Internal(&detM, M);
964: *vol = 0.5*detM;
965: (void)PetscLogFlops(5.0);
966: }
968: PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[])
969: {
970: DMPlex_Det2D_Internal(vol, coords);
971: *vol *= 0.5;
972: }
974: PETSC_UNUSED
975: PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[])
976: {
977: /* Signed volume is 1/6th of the determinant
979: | 1 1 1 1 |
980: | x0 x1 x2 x3 |
981: | y0 y1 y2 y3 |
982: | z0 z1 z2 z3 |
984: but if x0,y0,z0 is the origin, we have
986: | x1 x2 x3 |
987: | y1 y2 y3 |
988: | z1 z2 z3 |
989: */
990: const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2];
991: const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2];
992: const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2];
993: const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
994: PetscReal M[9], detM;
995: M[0] = x1; M[1] = x2; M[2] = x3;
996: M[3] = y1; M[4] = y2; M[5] = y3;
997: M[6] = z1; M[7] = z2; M[8] = z3;
998: DMPlex_Det3D_Internal(&detM, M);
999: *vol = -onesixth*detM;
1000: (void)PetscLogFlops(10.0);
1001: }
1003: PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[])
1004: {
1005: const PetscReal onesixth = ((PetscReal)1./(PetscReal)6.);
1006: DMPlex_Det3D_Internal(vol, coords);
1007: *vol *= -onesixth;
1008: }
1010: static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1011: {
1012: PetscSection coordSection;
1013: Vec coordinates;
1014: const PetscScalar *coords;
1015: PetscInt dim, d, off;
1019: DMGetCoordinatesLocal(dm, &coordinates);
1020: DMGetCoordinateSection(dm, &coordSection);
1021: PetscSectionGetDof(coordSection,e,&dim);
1022: if (!dim) return(0);
1023: PetscSectionGetOffset(coordSection,e,&off);
1024: VecGetArrayRead(coordinates,&coords);
1025: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);}
1026: VecRestoreArrayRead(coordinates,&coords);
1027: *detJ = 1.;
1028: if (J) {
1029: for (d = 0; d < dim * dim; d++) J[d] = 0.;
1030: for (d = 0; d < dim; d++) J[d * dim + d] = 1.;
1031: if (invJ) {
1032: for (d = 0; d < dim * dim; d++) invJ[d] = 0.;
1033: for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.;
1034: }
1035: }
1036: return(0);
1037: }
1039: static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1040: {
1041: PetscSection coordSection;
1042: Vec coordinates;
1043: PetscScalar *coords = NULL;
1044: PetscInt numCoords, d, pStart, pEnd, numSelfCoords = 0;
1048: DMGetCoordinatesLocal(dm, &coordinates);
1049: DMGetCoordinateSection(dm, &coordSection);
1050: PetscSectionGetChart(coordSection,&pStart,&pEnd);
1051: if (e >= pStart && e < pEnd) {PetscSectionGetDof(coordSection,e,&numSelfCoords);}
1052: DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1053: numCoords = numSelfCoords ? numSelfCoords : numCoords;
1054: if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL");
1055: *detJ = 0.0;
1056: if (numCoords == 6) {
1057: const PetscInt dim = 3;
1058: PetscReal R[9], J0;
1060: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1061: DMPlexComputeProjection3Dto1D(coords, R);
1062: if (J) {
1063: J0 = 0.5*PetscRealPart(coords[1]);
1064: J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2];
1065: J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5];
1066: J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8];
1067: DMPlex_Det3D_Internal(detJ, J);
1068: if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1069: }
1070: } else if (numCoords == 4) {
1071: const PetscInt dim = 2;
1072: PetscReal R[4], J0;
1074: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1075: DMPlexComputeProjection2Dto1D(coords, R);
1076: if (J) {
1077: J0 = 0.5*PetscRealPart(coords[1]);
1078: J[0] = R[0]*J0; J[1] = R[1];
1079: J[2] = R[2]*J0; J[3] = R[3];
1080: DMPlex_Det2D_Internal(detJ, J);
1081: if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1082: }
1083: } else if (numCoords == 2) {
1084: const PetscInt dim = 1;
1086: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1087: if (J) {
1088: J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0]));
1089: *detJ = J[0];
1090: PetscLogFlops(2.0);
1091: if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);}
1092: }
1093: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords);
1094: DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1095: return(0);
1096: }
1098: static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1099: {
1100: PetscSection coordSection;
1101: Vec coordinates;
1102: PetscScalar *coords = NULL;
1103: PetscInt numCoords, d, f, g;
1107: DMGetCoordinatesLocal(dm, &coordinates);
1108: DMGetCoordinateSection(dm, &coordSection);
1109: DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1110: *detJ = 0.0;
1111: if (numCoords == 9) {
1112: const PetscInt dim = 3;
1113: PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
1115: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1116: DMPlexComputeProjection3Dto2D(numCoords, coords, R);
1117: if (J) {
1118: const PetscInt pdim = 2;
1120: for (d = 0; d < pdim; d++) {
1121: for (f = 0; f < pdim; f++) {
1122: J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1123: }
1124: }
1125: PetscLogFlops(8.0);
1126: DMPlex_Det3D_Internal(detJ, J0);
1127: for (d = 0; d < dim; d++) {
1128: for (f = 0; f < dim; f++) {
1129: J[d*dim+f] = 0.0;
1130: for (g = 0; g < dim; g++) {
1131: J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
1132: }
1133: }
1134: }
1135: PetscLogFlops(18.0);
1136: }
1137: if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1138: } else if (numCoords == 6) {
1139: const PetscInt dim = 2;
1141: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1142: if (J) {
1143: for (d = 0; d < dim; d++) {
1144: for (f = 0; f < dim; f++) {
1145: J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d]));
1146: }
1147: }
1148: PetscLogFlops(8.0);
1149: DMPlex_Det2D_Internal(detJ, J);
1150: }
1151: if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1152: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords);
1153: DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1154: return(0);
1155: }
1157: static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1158: {
1159: PetscSection coordSection;
1160: Vec coordinates;
1161: PetscScalar *coords = NULL;
1162: PetscInt numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd;
1166: DMGetCoordinatesLocal(dm, &coordinates);
1167: DMGetCoordinateSection(dm, &coordSection);
1168: PetscSectionGetChart(coordSection,&pStart,&pEnd);
1169: if (e >= pStart && e < pEnd) {PetscSectionGetDof(coordSection,e,&numSelfCoords);}
1170: DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1171: numCoords = numSelfCoords ? numSelfCoords : numCoords;
1172: if (!Nq) {
1173: *detJ = 0.0;
1174: if (numCoords == 12) {
1175: const PetscInt dim = 3;
1176: PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0};
1178: if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1179: DMPlexComputeProjection3Dto2D(numCoords, coords, R);
1180: if (J) {
1181: const PetscInt pdim = 2;
1183: for (d = 0; d < pdim; d++) {
1184: J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1185: J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d]));
1186: }
1187: PetscLogFlops(8.0);
1188: DMPlex_Det3D_Internal(detJ, J0);
1189: for (d = 0; d < dim; d++) {
1190: for (f = 0; f < dim; f++) {
1191: J[d*dim+f] = 0.0;
1192: for (g = 0; g < dim; g++) {
1193: J[d*dim+f] += R[d*dim+g]*J0[g*dim+f];
1194: }
1195: }
1196: }
1197: PetscLogFlops(18.0);
1198: }
1199: if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1200: } else if (numCoords == 8) {
1201: const PetscInt dim = 2;
1203: if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1204: if (J) {
1205: for (d = 0; d < dim; d++) {
1206: J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1207: J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1208: }
1209: PetscLogFlops(8.0);
1210: DMPlex_Det2D_Internal(detJ, J);
1211: }
1212: if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);}
1213: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1214: } else {
1215: const PetscInt Nv = 4;
1216: const PetscInt dimR = 2;
1217: const PetscInt zToPlex[4] = {0, 1, 3, 2};
1218: PetscReal zOrder[12];
1219: PetscReal zCoeff[12];
1220: PetscInt i, j, k, l, dim;
1222: if (numCoords == 12) {
1223: dim = 3;
1224: } else if (numCoords == 8) {
1225: dim = 2;
1226: } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords);
1227: for (i = 0; i < Nv; i++) {
1228: PetscInt zi = zToPlex[i];
1230: for (j = 0; j < dim; j++) {
1231: zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1232: }
1233: }
1234: for (j = 0; j < dim; j++) {
1235: zCoeff[dim * 0 + j] = 0.25 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1236: zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1237: zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1238: zCoeff[dim * 3 + j] = 0.25 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]);
1239: }
1240: for (i = 0; i < Nq; i++) {
1241: PetscReal xi = points[dimR * i], eta = points[dimR * i + 1];
1243: if (v) {
1244: PetscReal extPoint[4];
1246: extPoint[0] = 1.;
1247: extPoint[1] = xi;
1248: extPoint[2] = eta;
1249: extPoint[3] = xi * eta;
1250: for (j = 0; j < dim; j++) {
1251: PetscReal val = 0.;
1253: for (k = 0; k < Nv; k++) {
1254: val += extPoint[k] * zCoeff[dim * k + j];
1255: }
1256: v[i * dim + j] = val;
1257: }
1258: }
1259: if (J) {
1260: PetscReal extJ[8];
1262: extJ[0] = 0.;
1263: extJ[1] = 0.;
1264: extJ[2] = 1.;
1265: extJ[3] = 0.;
1266: extJ[4] = 0.;
1267: extJ[5] = 1.;
1268: extJ[6] = eta;
1269: extJ[7] = xi;
1270: for (j = 0; j < dim; j++) {
1271: for (k = 0; k < dimR; k++) {
1272: PetscReal val = 0.;
1274: for (l = 0; l < Nv; l++) {
1275: val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1276: }
1277: J[i * dim * dim + dim * j + k] = val;
1278: }
1279: }
1280: if (dim == 3) { /* put the cross product in the third component of the Jacobian */
1281: PetscReal x, y, z;
1282: PetscReal *iJ = &J[i * dim * dim];
1283: PetscReal norm;
1285: x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0];
1286: y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1];
1287: z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0];
1288: norm = PetscSqrtReal(x * x + y * y + z * z);
1289: iJ[2] = x / norm;
1290: iJ[5] = y / norm;
1291: iJ[8] = z / norm;
1292: DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1293: if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1294: } else {
1295: DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]);
1296: if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1297: }
1298: }
1299: }
1300: }
1301: DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);
1302: return(0);
1303: }
1305: static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1306: {
1307: PetscSection coordSection;
1308: Vec coordinates;
1309: PetscScalar *coords = NULL;
1310: const PetscInt dim = 3;
1311: PetscInt d;
1315: DMGetCoordinatesLocal(dm, &coordinates);
1316: DMGetCoordinateSection(dm, &coordSection);
1317: DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);
1318: *detJ = 0.0;
1319: if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);}
1320: if (J) {
1321: for (d = 0; d < dim; d++) {
1322: /* I orient with outward face normals */
1323: J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d]));
1324: J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1325: J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1326: }
1327: PetscLogFlops(18.0);
1328: DMPlex_Det3D_Internal(detJ, J);
1329: }
1330: if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1331: DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);
1332: return(0);
1333: }
1335: static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1336: {
1337: PetscSection coordSection;
1338: Vec coordinates;
1339: PetscScalar *coords = NULL;
1340: const PetscInt dim = 3;
1341: PetscInt d;
1345: DMGetCoordinatesLocal(dm, &coordinates);
1346: DMGetCoordinateSection(dm, &coordSection);
1347: DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);
1348: if (!Nq) {
1349: *detJ = 0.0;
1350: if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);}
1351: if (J) {
1352: for (d = 0; d < dim; d++) {
1353: J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d]));
1354: J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d]));
1355: J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d]));
1356: }
1357: PetscLogFlops(18.0);
1358: DMPlex_Det3D_Internal(detJ, J);
1359: }
1360: if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);}
1361: } else {
1362: const PetscInt Nv = 8;
1363: const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
1364: const PetscInt dim = 3;
1365: const PetscInt dimR = 3;
1366: PetscReal zOrder[24];
1367: PetscReal zCoeff[24];
1368: PetscInt i, j, k, l;
1370: for (i = 0; i < Nv; i++) {
1371: PetscInt zi = zToPlex[i];
1373: for (j = 0; j < dim; j++) {
1374: zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]);
1375: }
1376: }
1377: for (j = 0; j < dim; j++) {
1378: zCoeff[dim * 0 + j] = 0.125 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1379: zCoeff[dim * 1 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1380: zCoeff[dim * 2 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1381: zCoeff[dim * 3 + j] = 0.125 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1382: zCoeff[dim * 4 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1383: zCoeff[dim * 5 + j] = 0.125 * (+ zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1384: zCoeff[dim * 6 + j] = 0.125 * (+ zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1385: zCoeff[dim * 7 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]);
1386: }
1387: for (i = 0; i < Nq; i++) {
1388: PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2];
1390: if (v) {
1391: PetscReal extPoint[8];
1393: extPoint[0] = 1.;
1394: extPoint[1] = xi;
1395: extPoint[2] = eta;
1396: extPoint[3] = xi * eta;
1397: extPoint[4] = theta;
1398: extPoint[5] = theta * xi;
1399: extPoint[6] = theta * eta;
1400: extPoint[7] = theta * eta * xi;
1401: for (j = 0; j < dim; j++) {
1402: PetscReal val = 0.;
1404: for (k = 0; k < Nv; k++) {
1405: val += extPoint[k] * zCoeff[dim * k + j];
1406: }
1407: v[i * dim + j] = val;
1408: }
1409: }
1410: if (J) {
1411: PetscReal extJ[24];
1413: extJ[0] = 0. ; extJ[1] = 0. ; extJ[2] = 0. ;
1414: extJ[3] = 1. ; extJ[4] = 0. ; extJ[5] = 0. ;
1415: extJ[6] = 0. ; extJ[7] = 1. ; extJ[8] = 0. ;
1416: extJ[9] = eta ; extJ[10] = xi ; extJ[11] = 0. ;
1417: extJ[12] = 0. ; extJ[13] = 0. ; extJ[14] = 1. ;
1418: extJ[15] = theta ; extJ[16] = 0. ; extJ[17] = xi ;
1419: extJ[18] = 0. ; extJ[19] = theta ; extJ[20] = eta ;
1420: extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi;
1422: for (j = 0; j < dim; j++) {
1423: for (k = 0; k < dimR; k++) {
1424: PetscReal val = 0.;
1426: for (l = 0; l < Nv; l++) {
1427: val += zCoeff[dim * l + j] * extJ[dimR * l + k];
1428: }
1429: J[i * dim * dim + dim * j + k] = val;
1430: }
1431: }
1432: DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]);
1433: if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);}
1434: }
1435: }
1436: }
1437: DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);
1438: return(0);
1439: }
1441: static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1442: {
1443: PetscInt depth, dim, coordDim, coneSize, i;
1444: PetscInt Nq = 0;
1445: const PetscReal *points = NULL;
1446: DMLabel depthLabel;
1447: PetscReal v0[3] = {-1.}, J0[9] = {-1.}, detJ0 = -1.;
1448: PetscBool isAffine = PETSC_TRUE;
1449: PetscErrorCode ierr;
1452: DMPlexGetDepth(dm, &depth);
1453: DMPlexGetConeSize(dm, cell, &coneSize);
1454: DMPlexGetDepthLabel(dm, &depthLabel);
1455: DMLabelGetValue(depthLabel, cell, &dim);
1456: if (depth == 1 && dim == 1) {
1457: DMGetDimension(dm, &dim);
1458: }
1459: DMGetCoordinateDim(dm, &coordDim);
1460: if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim);
1461: if (quad) {PetscQuadratureGetData(quad, NULL, NULL, &Nq, &points, NULL);}
1462: switch (dim) {
1463: case 0:
1464: DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);
1465: isAffine = PETSC_FALSE;
1466: break;
1467: case 1:
1468: if (Nq) {
1469: DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);
1470: } else {
1471: DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);
1472: }
1473: break;
1474: case 2:
1475: switch (coneSize) {
1476: case 3:
1477: if (Nq) {
1478: DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);
1479: } else {
1480: DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);
1481: }
1482: break;
1483: case 4:
1484: DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);
1485: isAffine = PETSC_FALSE;
1486: break;
1487: default:
1488: SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1489: }
1490: break;
1491: case 3:
1492: switch (coneSize) {
1493: case 4:
1494: if (Nq) {
1495: DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);
1496: } else {
1497: DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);
1498: }
1499: break;
1500: case 6: /* Faces */
1501: case 8: /* Vertices */
1502: DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);
1503: isAffine = PETSC_FALSE;
1504: break;
1505: default:
1506: SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell);
1507: }
1508: break;
1509: default:
1510: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1511: }
1512: if (isAffine && Nq) {
1513: if (v) {
1514: for (i = 0; i < Nq; i++) {
1515: CoordinatesRefToReal(coordDim,dim,v0, J0, &points[dim * i], &v[coordDim * i]);
1516: }
1517: }
1518: if (detJ) {
1519: for (i = 0; i < Nq; i++) {
1520: detJ[i] = detJ0;
1521: }
1522: }
1523: if (J) {
1524: PetscInt k;
1526: for (i = 0, k = 0; i < Nq; i++) {
1527: PetscInt j;
1529: for (j = 0; j < coordDim * coordDim; j++, k++) {
1530: J[k] = J0[j];
1531: }
1532: }
1533: }
1534: if (invJ) {
1535: PetscInt k;
1536: switch (coordDim) {
1537: case 0:
1538: break;
1539: case 1:
1540: invJ[0] = 1./J0[0];
1541: break;
1542: case 2:
1543: DMPlex_Invert2D_Internal(invJ, J0, detJ0);
1544: break;
1545: case 3:
1546: DMPlex_Invert3D_Internal(invJ, J0, detJ0);
1547: break;
1548: }
1549: for (i = 1, k = coordDim * coordDim; i < Nq; i++) {
1550: PetscInt j;
1552: for (j = 0; j < coordDim * coordDim; j++, k++) {
1553: invJ[k] = invJ[j];
1554: }
1555: }
1556: }
1557: }
1558: return(0);
1559: }
1561: /*@C
1562: DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell
1564: Collective on DM
1566: Input Arguments:
1567: + dm - the DM
1568: - cell - the cell
1570: Output Arguments:
1571: + v0 - the translation part of this affine transform
1572: . J - the Jacobian of the transform from the reference element
1573: . invJ - the inverse of the Jacobian
1574: - detJ - the Jacobian determinant
1576: Level: advanced
1578: Fortran Notes:
1579: Since it returns arrays, this routine is only available in Fortran 90, and you must
1580: include petsc.h90 in your code.
1582: .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec()
1583: @*/
1584: PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1585: {
1589: DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);
1590: return(0);
1591: }
1593: static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
1594: {
1595: PetscQuadrature feQuad;
1596: PetscSection coordSection;
1597: Vec coordinates;
1598: PetscScalar *coords = NULL;
1599: const PetscReal *quadPoints;
1600: PetscReal *basisDer, *basis, detJt;
1601: PetscInt dim, cdim, pdim, qdim, Nq, numCoords, q;
1602: PetscErrorCode ierr;
1605: DMGetCoordinatesLocal(dm, &coordinates);
1606: DMGetCoordinateSection(dm, &coordSection);
1607: DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);
1608: DMGetDimension(dm, &dim);
1609: DMGetCoordinateDim(dm, &cdim);
1610: if (!quad) { /* use the first point of the first functional of the dual space */
1611: PetscDualSpace dsp;
1613: PetscFEGetDualSpace(fe, &dsp);
1614: PetscDualSpaceGetFunctional(dsp, 0, &quad);
1615: PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);
1616: Nq = 1;
1617: } else {
1618: PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL);
1619: }
1620: PetscFEGetDimension(fe, &pdim);
1621: PetscFEGetQuadrature(fe, &feQuad);
1622: if (feQuad == quad) {
1623: PetscFEGetDefaultTabulation(fe, &basis, J ? &basisDer : NULL, NULL);
1624: if (numCoords != pdim*cdim) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "There are %d coordinates for point %d != %d*%d", numCoords, point, pdim, cdim);
1625: } else {
1626: PetscFEGetTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);
1627: }
1628: if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim);
1629: if (v) {
1630: PetscMemzero(v, Nq*cdim*sizeof(PetscReal));
1631: for (q = 0; q < Nq; ++q) {
1632: PetscInt i, k;
1634: for (k = 0; k < pdim; ++k)
1635: for (i = 0; i < cdim; ++i)
1636: v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]);
1637: PetscLogFlops(2.0*pdim*cdim);
1638: }
1639: }
1640: if (J) {
1641: PetscMemzero(J, Nq*cdim*cdim*sizeof(PetscReal));
1642: for (q = 0; q < Nq; ++q) {
1643: PetscInt i, j, k, c, r;
1645: /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */
1646: for (k = 0; k < pdim; ++k)
1647: for (j = 0; j < dim; ++j)
1648: for (i = 0; i < cdim; ++i)
1649: J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]);
1650: PetscLogFlops(2.0*pdim*dim*cdim);
1651: if (cdim > dim) {
1652: for (c = dim; c < cdim; ++c)
1653: for (r = 0; r < cdim; ++r)
1654: J[r*cdim+c] = r == c ? 1.0 : 0.0;
1655: }
1656: if (!detJ && !invJ) continue;
1657: detJt = 0.;
1658: switch (cdim) {
1659: case 3:
1660: DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]);
1661: if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
1662: break;
1663: case 2:
1664: DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]);
1665: if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);}
1666: break;
1667: case 1:
1668: detJt = J[q*cdim*dim];
1669: if (invJ) invJ[q*cdim*dim] = 1.0/detJt;
1670: }
1671: if (detJ) detJ[q] = detJt;
1672: }
1673: }
1674: else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ");
1675: if (feQuad != quad) {
1676: PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, J ? &basisDer : NULL, NULL);
1677: }
1678: DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);
1679: return(0);
1680: }
1682: /*@C
1683: DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell
1685: Collective on DM
1687: Input Arguments:
1688: + dm - the DM
1689: . cell - the cell
1690: - quad - the quadrature containing the points in the reference element where the geometry will be evaluated. If quad == NULL, geometry will be
1691: evaluated at the first vertex of the reference element
1693: Output Arguments:
1694: + v - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element
1695: . J - the Jacobian of the transform from the reference element at each quadrature point
1696: . invJ - the inverse of the Jacobian at each quadrature point
1697: - detJ - the Jacobian determinant at each quadrature point
1699: Level: advanced
1701: Fortran Notes:
1702: Since it returns arrays, this routine is only available in Fortran 90, and you must
1703: include petsc.h90 in your code.
1705: .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
1706: @*/
1707: PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ)
1708: {
1709: PetscFE fe = NULL;
1713: if (dm->coordinateDM) {
1714: PetscClassId id;
1715: PetscInt numFields;
1716: PetscDS prob = dm->coordinateDM->prob;
1717: PetscObject disc;
1719: PetscDSGetNumFields(prob, &numFields);
1720: if (numFields) {
1721: PetscDSGetDiscretization(prob,0,&disc);
1722: PetscObjectGetClassId(disc,&id);
1723: if (id == PETSCFE_CLASSID) {
1724: fe = (PetscFE) disc;
1725: }
1726: }
1727: }
1728: if (!fe) {DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);}
1729: else {DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);}
1730: return(0);
1731: }
1733: static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1734: {
1735: PetscSection coordSection;
1736: Vec coordinates;
1737: PetscScalar *coords = NULL;
1738: PetscScalar tmp[2];
1739: PetscInt coordSize;
1743: DMGetCoordinatesLocal(dm, &coordinates);
1744: DMGetCoordinateSection(dm, &coordSection);
1745: DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);
1746: if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now");
1747: DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);
1748: if (centroid) {
1749: centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]);
1750: centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]);
1751: }
1752: if (normal) {
1753: PetscReal norm;
1755: normal[0] = -PetscRealPart(coords[1] - tmp[1]);
1756: normal[1] = PetscRealPart(coords[0] - tmp[0]);
1757: norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]);
1758: normal[0] /= norm;
1759: normal[1] /= norm;
1760: }
1761: if (vol) {
1762: *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1])));
1763: }
1764: DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);
1765: return(0);
1766: }
1768: /* Centroid_i = (\sum_n A_n Cn_i ) / A */
1769: static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1770: {
1771: PetscSection coordSection;
1772: Vec coordinates;
1773: PetscScalar *coords = NULL;
1774: PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9];
1775: PetscInt tdim = 2, coordSize, numCorners, p, d, e;
1779: DMGetCoordinatesLocal(dm, &coordinates);
1780: DMPlexGetConeSize(dm, cell, &numCorners);
1781: DMGetCoordinateSection(dm, &coordSection);
1782: DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);
1783: DMGetCoordinateDim(dm, &dim);
1784: if (dim > 2 && centroid) {
1785: v0[0] = PetscRealPart(coords[0]);
1786: v0[1] = PetscRealPart(coords[1]);
1787: v0[2] = PetscRealPart(coords[2]);
1788: }
1789: if (normal) {
1790: if (dim > 2) {
1791: const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]);
1792: const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]);
1793: const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]);
1794: PetscReal norm;
1796: normal[0] = y0*z1 - z0*y1;
1797: normal[1] = z0*x1 - x0*z1;
1798: normal[2] = x0*y1 - y0*x1;
1799: norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]);
1800: normal[0] /= norm;
1801: normal[1] /= norm;
1802: normal[2] /= norm;
1803: } else {
1804: for (d = 0; d < dim; ++d) normal[d] = 0.0;
1805: }
1806: }
1807: if (dim == 3) {DMPlexComputeProjection3Dto2D(coordSize, coords, R);}
1808: for (p = 0; p < numCorners; ++p) {
1809: /* Need to do this copy to get types right */
1810: for (d = 0; d < tdim; ++d) {
1811: ctmp[d] = PetscRealPart(coords[p*tdim+d]);
1812: ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]);
1813: }
1814: Volume_Triangle_Origin_Internal(&vtmp, ctmp);
1815: vsum += vtmp;
1816: for (d = 0; d < tdim; ++d) {
1817: csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp;
1818: }
1819: }
1820: for (d = 0; d < tdim; ++d) {
1821: csum[d] /= (tdim+1)*vsum;
1822: }
1823: DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);
1824: if (vol) *vol = PetscAbsReal(vsum);
1825: if (centroid) {
1826: if (dim > 2) {
1827: for (d = 0; d < dim; ++d) {
1828: centroid[d] = v0[d];
1829: for (e = 0; e < dim; ++e) {
1830: centroid[d] += R[d*dim+e]*csum[e];
1831: }
1832: }
1833: } else for (d = 0; d < dim; ++d) centroid[d] = csum[d];
1834: }
1835: return(0);
1836: }
1838: /* Centroid_i = (\sum_n V_n Cn_i ) / V */
1839: static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1840: {
1841: PetscSection coordSection;
1842: Vec coordinates;
1843: PetscScalar *coords = NULL;
1844: PetscReal vsum = 0.0, vtmp, coordsTmp[3*3];
1845: const PetscInt *faces, *facesO;
1846: PetscInt numFaces, f, coordSize, numCorners, p, d;
1847: PetscErrorCode ierr;
1850: if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim);
1851: DMGetCoordinatesLocal(dm, &coordinates);
1852: DMGetCoordinateSection(dm, &coordSection);
1854: if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0;
1855: DMPlexGetConeSize(dm, cell, &numFaces);
1856: DMPlexGetCone(dm, cell, &faces);
1857: DMPlexGetConeOrientation(dm, cell, &facesO);
1858: for (f = 0; f < numFaces; ++f) {
1859: DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);
1860: numCorners = coordSize/dim;
1861: switch (numCorners) {
1862: case 3:
1863: for (d = 0; d < dim; ++d) {
1864: coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
1865: coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
1866: coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]);
1867: }
1868: Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1869: if (facesO[f] < 0) vtmp = -vtmp;
1870: vsum += vtmp;
1871: if (centroid) { /* Centroid of OABC = (a+b+c)/4 */
1872: for (d = 0; d < dim; ++d) {
1873: for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1874: }
1875: }
1876: break;
1877: case 4:
1878: /* DO FOR PYRAMID */
1879: /* First tet */
1880: for (d = 0; d < dim; ++d) {
1881: coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]);
1882: coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]);
1883: coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
1884: }
1885: Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1886: if (facesO[f] < 0) vtmp = -vtmp;
1887: vsum += vtmp;
1888: if (centroid) {
1889: for (d = 0; d < dim; ++d) {
1890: for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1891: }
1892: }
1893: /* Second tet */
1894: for (d = 0; d < dim; ++d) {
1895: coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]);
1896: coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]);
1897: coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]);
1898: }
1899: Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp);
1900: if (facesO[f] < 0) vtmp = -vtmp;
1901: vsum += vtmp;
1902: if (centroid) {
1903: for (d = 0; d < dim; ++d) {
1904: for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp;
1905: }
1906: }
1907: break;
1908: default:
1909: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners);
1910: }
1911: DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);
1912: }
1913: if (vol) *vol = PetscAbsReal(vsum);
1914: if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0;
1915: if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4);
1916: return(0);
1917: }
1919: /*@C
1920: DMPlexComputeCellGeometryFVM - Compute the volume for a given cell
1922: Collective on DM
1924: Input Arguments:
1925: + dm - the DM
1926: - cell - the cell
1928: Output Arguments:
1929: + volume - the cell volume
1930: . centroid - the cell centroid
1931: - normal - the cell normal, if appropriate
1933: Level: advanced
1935: Fortran Notes:
1936: Since it returns arrays, this routine is only available in Fortran 90, and you must
1937: include petsc.h90 in your code.
1939: .seealso: DMGetCoordinateSection(), DMGetCoordinateVec()
1940: @*/
1941: PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[])
1942: {
1943: PetscInt depth, dim;
1947: DMPlexGetDepth(dm, &depth);
1948: DMGetDimension(dm, &dim);
1949: if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
1950: /* We need to keep a pointer to the depth label */
1951: DMGetLabelValue(dm, "depth", cell, &depth);
1952: /* Cone size is now the number of faces */
1953: switch (depth) {
1954: case 1:
1955: DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);
1956: break;
1957: case 2:
1958: DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);
1959: break;
1960: case 3:
1961: DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);
1962: break;
1963: default:
1964: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim);
1965: }
1966: return(0);
1967: }
1969: /*@
1970: DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh
1972: Collective on dm
1974: Input Parameter:
1975: . dm - The DMPlex
1977: Output Parameter:
1978: . cellgeom - A vector with the cell geometry data for each cell
1980: Level: beginner
1982: .keywords: DMPlexComputeCellGeometryFEM()
1983: @*/
1984: PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom)
1985: {
1986: DM dmCell;
1987: Vec coordinates;
1988: PetscSection coordSection, sectionCell;
1989: PetscScalar *cgeom;
1990: PetscInt cStart, cEnd, cMax, c;
1994: DMClone(dm, &dmCell);
1995: DMGetCoordinateSection(dm, &coordSection);
1996: DMGetCoordinatesLocal(dm, &coordinates);
1997: DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);
1998: DMSetCoordinatesLocal(dmCell, coordinates);
1999: PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);
2000: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
2001: DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
2002: cEnd = cMax < 0 ? cEnd : cMax;
2003: PetscSectionSetChart(sectionCell, cStart, cEnd);
2004: /* TODO This needs to be multiplied by Nq for non-affine */
2005: for (c = cStart; c < cEnd; ++c) {PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));}
2006: PetscSectionSetUp(sectionCell);
2007: DMSetDefaultSection(dmCell, sectionCell);
2008: PetscSectionDestroy(§ionCell);
2009: DMCreateLocalVector(dmCell, cellgeom);
2010: VecGetArray(*cellgeom, &cgeom);
2011: for (c = cStart; c < cEnd; ++c) {
2012: PetscFECellGeom *cg;
2014: DMPlexPointLocalRef(dmCell, c, cgeom, &cg);
2015: PetscMemzero(cg, sizeof(*cg));
2016: DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);
2017: if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c);
2018: }
2019: VecRestoreArray(*cellgeom, &cgeom);
2020: DMDestroy(&dmCell);
2021: return(0);
2022: }
2024: /*@
2025: DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method
2027: Input Parameter:
2028: . dm - The DM
2030: Output Parameters:
2031: + cellgeom - A Vec of PetscFVCellGeom data
2032: . facegeom - A Vec of PetscFVFaceGeom data
2034: Level: developer
2036: .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM()
2037: @*/
2038: PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom)
2039: {
2040: DM dmFace, dmCell;
2041: DMLabel ghostLabel;
2042: PetscSection sectionFace, sectionCell;
2043: PetscSection coordSection;
2044: Vec coordinates;
2045: PetscScalar *fgeom, *cgeom;
2046: PetscReal minradius, gminradius;
2047: PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
2051: DMGetDimension(dm, &dim);
2052: DMGetCoordinateSection(dm, &coordSection);
2053: DMGetCoordinatesLocal(dm, &coordinates);
2054: /* Make cell centroids and volumes */
2055: DMClone(dm, &dmCell);
2056: DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);
2057: DMSetCoordinatesLocal(dmCell, coordinates);
2058: PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);
2059: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
2060: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
2061: PetscSectionSetChart(sectionCell, cStart, cEnd);
2062: for (c = cStart; c < cEnd; ++c) {PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));}
2063: PetscSectionSetUp(sectionCell);
2064: DMSetDefaultSection(dmCell, sectionCell);
2065: PetscSectionDestroy(§ionCell);
2066: DMCreateLocalVector(dmCell, cellgeom);
2067: if (cEndInterior < 0) {
2068: cEndInterior = cEnd;
2069: }
2070: VecGetArray(*cellgeom, &cgeom);
2071: for (c = cStart; c < cEndInterior; ++c) {
2072: PetscFVCellGeom *cg;
2074: DMPlexPointLocalRef(dmCell, c, cgeom, &cg);
2075: PetscMemzero(cg, sizeof(*cg));
2076: DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);
2077: }
2078: /* Compute face normals and minimum cell radius */
2079: DMClone(dm, &dmFace);
2080: PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);
2081: DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
2082: PetscSectionSetChart(sectionFace, fStart, fEnd);
2083: for (f = fStart; f < fEnd; ++f) {PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));}
2084: PetscSectionSetUp(sectionFace);
2085: DMSetDefaultSection(dmFace, sectionFace);
2086: PetscSectionDestroy(§ionFace);
2087: DMCreateLocalVector(dmFace, facegeom);
2088: VecGetArray(*facegeom, &fgeom);
2089: DMGetLabel(dm, "ghost", &ghostLabel);
2090: minradius = PETSC_MAX_REAL;
2091: for (f = fStart; f < fEnd; ++f) {
2092: PetscFVFaceGeom *fg;
2093: PetscReal area;
2094: PetscInt ghost = -1, d, numChildren;
2096: if (ghostLabel) {DMLabelGetValue(ghostLabel, f, &ghost);}
2097: DMPlexGetTreeChildren(dm,f,&numChildren,NULL);
2098: if (ghost >= 0 || numChildren) continue;
2099: DMPlexPointLocalRef(dmFace, f, fgeom, &fg);
2100: DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);
2101: for (d = 0; d < dim; ++d) fg->normal[d] *= area;
2102: /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
2103: {
2104: PetscFVCellGeom *cL, *cR;
2105: PetscInt ncells;
2106: const PetscInt *cells;
2107: PetscReal *lcentroid, *rcentroid;
2108: PetscReal l[3], r[3], v[3];
2110: DMPlexGetSupport(dm, f, &cells);
2111: DMPlexGetSupportSize(dm, f, &ncells);
2112: DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);
2113: lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
2114: if (ncells > 1) {
2115: DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);
2116: rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
2117: }
2118: else {
2119: rcentroid = fg->centroid;
2120: }
2121: DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);
2122: DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);
2123: DMPlex_WaxpyD_Internal(dim, -1, l, r, v);
2124: if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) {
2125: for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
2126: }
2127: if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) {
2128: if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]);
2129: if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]);
2130: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
2131: }
2132: if (cells[0] < cEndInterior) {
2133: DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v);
2134: minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2135: }
2136: if (ncells > 1 && cells[1] < cEndInterior) {
2137: DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v);
2138: minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v));
2139: }
2140: }
2141: }
2142: MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));
2143: DMPlexSetMinRadius(dm, gminradius);
2144: /* Compute centroids of ghost cells */
2145: for (c = cEndInterior; c < cEnd; ++c) {
2146: PetscFVFaceGeom *fg;
2147: const PetscInt *cone, *support;
2148: PetscInt coneSize, supportSize, s;
2150: DMPlexGetConeSize(dmCell, c, &coneSize);
2151: if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
2152: DMPlexGetCone(dmCell, c, &cone);
2153: DMPlexGetSupportSize(dmCell, cone[0], &supportSize);
2154: if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize);
2155: DMPlexGetSupport(dmCell, cone[0], &support);
2156: DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);
2157: for (s = 0; s < 2; ++s) {
2158: /* Reflect ghost centroid across plane of face */
2159: if (support[s] == c) {
2160: PetscFVCellGeom *ci;
2161: PetscFVCellGeom *cg;
2162: PetscReal c2f[3], a;
2164: DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);
2165: DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2166: a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal);
2167: DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);
2168: DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
2169: cg->volume = ci->volume;
2170: }
2171: }
2172: }
2173: VecRestoreArray(*facegeom, &fgeom);
2174: VecRestoreArray(*cellgeom, &cgeom);
2175: DMDestroy(&dmCell);
2176: DMDestroy(&dmFace);
2177: return(0);
2178: }
2180: /*@C
2181: DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face
2183: Not collective
2185: Input Argument:
2186: . dm - the DM
2188: Output Argument:
2189: . minradius - the minium cell radius
2191: Level: developer
2193: .seealso: DMGetCoordinates()
2194: @*/
2195: PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius)
2196: {
2200: *minradius = ((DM_Plex*) dm->data)->minradius;
2201: return(0);
2202: }
2204: /*@C
2205: DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face
2207: Logically collective
2209: Input Arguments:
2210: + dm - the DM
2211: - minradius - the minium cell radius
2213: Level: developer
2215: .seealso: DMSetCoordinates()
2216: @*/
2217: PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius)
2218: {
2221: ((DM_Plex*) dm->data)->minradius = minradius;
2222: return(0);
2223: }
2225: static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2226: {
2227: DMLabel ghostLabel;
2228: PetscScalar *dx, *grad, **gref;
2229: PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
2233: DMGetDimension(dm, &dim);
2234: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
2235: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
2236: DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);
2237: PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);
2238: DMGetLabel(dm, "ghost", &ghostLabel);
2239: PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);
2240: for (c = cStart; c < cEndInterior; c++) {
2241: const PetscInt *faces;
2242: PetscInt numFaces, usedFaces, f, d;
2243: PetscFVCellGeom *cg;
2244: PetscBool boundary;
2245: PetscInt ghost;
2247: DMPlexPointLocalRead(dmCell, c, cgeom, &cg);
2248: DMPlexGetConeSize(dm, c, &numFaces);
2249: DMPlexGetCone(dm, c, &faces);
2250: if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2251: for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2252: PetscFVCellGeom *cg1;
2253: PetscFVFaceGeom *fg;
2254: const PetscInt *fcells;
2255: PetscInt ncell, side;
2257: DMLabelGetValue(ghostLabel, faces[f], &ghost);
2258: DMIsBoundaryPoint(dm, faces[f], &boundary);
2259: if ((ghost >= 0) || boundary) continue;
2260: DMPlexGetSupport(dm, faces[f], &fcells);
2261: side = (c != fcells[0]); /* c is on left=0 or right=1 of face */
2262: ncell = fcells[!side]; /* the neighbor */
2263: DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);
2264: DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);
2265: for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
2266: gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */
2267: }
2268: if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
2269: PetscFVComputeGradient(fvm, usedFaces, dx, grad);
2270: for (f = 0, usedFaces = 0; f < numFaces; ++f) {
2271: DMLabelGetValue(ghostLabel, faces[f], &ghost);
2272: DMIsBoundaryPoint(dm, faces[f], &boundary);
2273: if ((ghost >= 0) || boundary) continue;
2274: for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
2275: ++usedFaces;
2276: }
2277: }
2278: PetscFree3(dx, grad, gref);
2279: return(0);
2280: }
2282: static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2283: {
2284: DMLabel ghostLabel;
2285: PetscScalar *dx, *grad, **gref;
2286: PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0;
2287: PetscSection neighSec;
2288: PetscInt (*neighbors)[2];
2289: PetscInt *counter;
2293: DMGetDimension(dm, &dim);
2294: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
2295: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
2296: if (cEndInterior < 0) {
2297: cEndInterior = cEnd;
2298: }
2299: PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);
2300: PetscSectionSetChart(neighSec,cStart,cEndInterior);
2301: DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
2302: DMGetLabel(dm, "ghost", &ghostLabel);
2303: for (f = fStart; f < fEnd; f++) {
2304: const PetscInt *fcells;
2305: PetscBool boundary;
2306: PetscInt ghost = -1;
2307: PetscInt numChildren, numCells, c;
2309: if (ghostLabel) {DMLabelGetValue(ghostLabel, f, &ghost);}
2310: DMIsBoundaryPoint(dm, f, &boundary);
2311: DMPlexGetTreeChildren(dm, f, &numChildren, NULL);
2312: if ((ghost >= 0) || boundary || numChildren) continue;
2313: DMPlexGetSupportSize(dm, f, &numCells);
2314: if (numCells == 2) {
2315: DMPlexGetSupport(dm, f, &fcells);
2316: for (c = 0; c < 2; c++) {
2317: PetscInt cell = fcells[c];
2319: if (cell >= cStart && cell < cEndInterior) {
2320: PetscSectionAddDof(neighSec,cell,1);
2321: }
2322: }
2323: }
2324: }
2325: PetscSectionSetUp(neighSec);
2326: PetscSectionGetMaxDof(neighSec,&maxNumFaces);
2327: PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);
2328: nStart = 0;
2329: PetscSectionGetStorageSize(neighSec,&nEnd);
2330: PetscMalloc1((nEnd-nStart),&neighbors);
2331: PetscCalloc1((cEndInterior-cStart),&counter);
2332: for (f = fStart; f < fEnd; f++) {
2333: const PetscInt *fcells;
2334: PetscBool boundary;
2335: PetscInt ghost = -1;
2336: PetscInt numChildren, numCells, c;
2338: if (ghostLabel) {DMLabelGetValue(ghostLabel, f, &ghost);}
2339: DMIsBoundaryPoint(dm, f, &boundary);
2340: DMPlexGetTreeChildren(dm, f, &numChildren, NULL);
2341: if ((ghost >= 0) || boundary || numChildren) continue;
2342: DMPlexGetSupportSize(dm, f, &numCells);
2343: if (numCells == 2) {
2344: DMPlexGetSupport(dm, f, &fcells);
2345: for (c = 0; c < 2; c++) {
2346: PetscInt cell = fcells[c], off;
2348: if (cell >= cStart && cell < cEndInterior) {
2349: PetscSectionGetOffset(neighSec,cell,&off);
2350: off += counter[cell - cStart]++;
2351: neighbors[off][0] = f;
2352: neighbors[off][1] = fcells[1 - c];
2353: }
2354: }
2355: }
2356: }
2357: PetscFree(counter);
2358: PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);
2359: for (c = cStart; c < cEndInterior; c++) {
2360: PetscInt numFaces, f, d, off, ghost = -1;
2361: PetscFVCellGeom *cg;
2363: DMPlexPointLocalRead(dmCell, c, cgeom, &cg);
2364: PetscSectionGetDof(neighSec, c, &numFaces);
2365: PetscSectionGetOffset(neighSec, c, &off);
2366: if (ghostLabel) {DMLabelGetValue(ghostLabel, c, &ghost);}
2367: if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces);
2368: for (f = 0; f < numFaces; ++f) {
2369: PetscFVCellGeom *cg1;
2370: PetscFVFaceGeom *fg;
2371: const PetscInt *fcells;
2372: PetscInt ncell, side, nface;
2374: nface = neighbors[off + f][0];
2375: ncell = neighbors[off + f][1];
2376: DMPlexGetSupport(dm,nface,&fcells);
2377: side = (c != fcells[0]);
2378: DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);
2379: DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);
2380: for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d];
2381: gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */
2382: }
2383: PetscFVComputeGradient(fvm, numFaces, dx, grad);
2384: for (f = 0; f < numFaces; ++f) {
2385: for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d];
2386: }
2387: }
2388: PetscFree3(dx, grad, gref);
2389: PetscSectionDestroy(&neighSec);
2390: PetscFree(neighbors);
2391: return(0);
2392: }
2394: /*@
2395: DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data
2397: Collective on DM
2399: Input Arguments:
2400: + dm - The DM
2401: . fvm - The PetscFV
2402: . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM()
2403: - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM()
2405: Output Parameters:
2406: + faceGeometry - The geometric factors for gradient calculation are inserted
2407: - dmGrad - The DM describing the layout of gradient data
2409: Level: developer
2411: .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM()
2412: @*/
2413: PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad)
2414: {
2415: DM dmFace, dmCell;
2416: PetscScalar *fgeom, *cgeom;
2417: PetscSection sectionGrad, parentSection;
2418: PetscInt dim, pdim, cStart, cEnd, cEndInterior, c;
2422: DMGetDimension(dm, &dim);
2423: PetscFVGetNumComponents(fvm, &pdim);
2424: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
2425: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
2426: /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
2427: VecGetDM(faceGeometry, &dmFace);
2428: VecGetDM(cellGeometry, &dmCell);
2429: VecGetArray(faceGeometry, &fgeom);
2430: VecGetArray(cellGeometry, &cgeom);
2431: DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);
2432: if (!parentSection) {
2433: BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);
2434: } else {
2435: BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);
2436: }
2437: VecRestoreArray(faceGeometry, &fgeom);
2438: VecRestoreArray(cellGeometry, &cgeom);
2439: /* Create storage for gradients */
2440: DMClone(dm, dmGrad);
2441: PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);
2442: PetscSectionSetChart(sectionGrad, cStart, cEnd);
2443: for (c = cStart; c < cEnd; ++c) {PetscSectionSetDof(sectionGrad, c, pdim*dim);}
2444: PetscSectionSetUp(sectionGrad);
2445: DMSetDefaultSection(*dmGrad, sectionGrad);
2446: PetscSectionDestroy(§ionGrad);
2447: return(0);
2448: }
2450: /*@
2451: DMPlexGetDataFVM - Retrieve precomputed cell geometry
2453: Collective on DM
2455: Input Arguments:
2456: + dm - The DM
2457: - fvm - The PetscFV
2459: Output Parameters:
2460: + cellGeometry - The cell geometry
2461: . faceGeometry - The face geometry
2462: - dmGrad - The gradient matrices
2464: Level: developer
2466: .seealso: DMPlexComputeGeometryFVM()
2467: @*/
2468: PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM)
2469: {
2470: PetscObject cellgeomobj, facegeomobj;
2474: PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);
2475: if (!cellgeomobj) {
2476: Vec cellgeomInt, facegeomInt;
2478: DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);
2479: PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);
2480: PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);
2481: VecDestroy(&cellgeomInt);
2482: VecDestroy(&facegeomInt);
2483: PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);
2484: }
2485: PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);
2486: if (cellgeom) *cellgeom = (Vec) cellgeomobj;
2487: if (facegeom) *facegeom = (Vec) facegeomobj;
2488: if (gradDM) {
2489: PetscObject gradobj;
2490: PetscBool computeGradients;
2492: PetscFVGetComputeGradients(fv,&computeGradients);
2493: if (!computeGradients) {
2494: *gradDM = NULL;
2495: return(0);
2496: }
2497: PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);
2498: if (!gradobj) {
2499: DM dmGradInt;
2501: DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);
2502: PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);
2503: DMDestroy(&dmGradInt);
2504: PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);
2505: }
2506: *gradDM = (DM) gradobj;
2507: }
2508: return(0);
2509: }
2511: static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work, PetscReal *resNeg, PetscReal *guess)
2512: {
2513: PetscInt l, m;
2516: if (dimC == dimR && dimR <= 3) {
2517: /* invert Jacobian, multiply */
2518: PetscScalar det, idet;
2520: switch (dimR) {
2521: case 1:
2522: invJ[0] = 1./ J[0];
2523: break;
2524: case 2:
2525: det = J[0] * J[3] - J[1] * J[2];
2526: idet = 1./det;
2527: invJ[0] = J[3] * idet;
2528: invJ[1] = -J[1] * idet;
2529: invJ[2] = -J[2] * idet;
2530: invJ[3] = J[0] * idet;
2531: break;
2532: case 3:
2533: {
2534: invJ[0] = J[4] * J[8] - J[5] * J[7];
2535: invJ[1] = J[2] * J[7] - J[1] * J[8];
2536: invJ[2] = J[1] * J[5] - J[2] * J[4];
2537: det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6];
2538: idet = 1./det;
2539: invJ[0] *= idet;
2540: invJ[1] *= idet;
2541: invJ[2] *= idet;
2542: invJ[3] = idet * (J[5] * J[6] - J[3] * J[8]);
2543: invJ[4] = idet * (J[0] * J[8] - J[2] * J[6]);
2544: invJ[5] = idet * (J[2] * J[3] - J[0] * J[5]);
2545: invJ[6] = idet * (J[3] * J[7] - J[4] * J[6]);
2546: invJ[7] = idet * (J[1] * J[6] - J[0] * J[7]);
2547: invJ[8] = idet * (J[0] * J[4] - J[1] * J[3]);
2548: }
2549: break;
2550: }
2551: for (l = 0; l < dimR; l++) {
2552: for (m = 0; m < dimC; m++) {
2553: guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m];
2554: }
2555: }
2556: } else {
2557: #if defined(PETSC_USE_COMPLEX)
2558: char transpose = 'C';
2559: #else
2560: char transpose = 'T';
2561: #endif
2562: PetscBLASInt m = dimR;
2563: PetscBLASInt n = dimC;
2564: PetscBLASInt one = 1;
2565: PetscBLASInt worksize = dimR * dimC, info;
2567: for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];}
2569: PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info));
2570: if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
2572: for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);}
2573: }
2574: return(0);
2575: }
2577: static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
2578: {
2579: PetscInt coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR);
2580: PetscScalar *coordsScalar = NULL;
2581: PetscReal *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg;
2582: PetscScalar *J, *invJ, *work;
2587: DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);
2588: if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
2589: DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);
2590: DMGetWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);
2591: cellCoords = &cellData[0];
2592: cellCoeffs = &cellData[coordSize];
2593: extJ = &cellData[2 * coordSize];
2594: resNeg = &cellData[2 * coordSize + dimR];
2595: invJ = &J[dimR * dimC];
2596: work = &J[2 * dimR * dimC];
2597: if (dimR == 2) {
2598: const PetscInt zToPlex[4] = {0, 1, 3, 2};
2600: for (i = 0; i < 4; i++) {
2601: PetscInt plexI = zToPlex[i];
2603: for (j = 0; j < dimC; j++) {
2604: cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2605: }
2606: }
2607: } else if (dimR == 3) {
2608: const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
2610: for (i = 0; i < 8; i++) {
2611: PetscInt plexI = zToPlex[i];
2613: for (j = 0; j < dimC; j++) {
2614: cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2615: }
2616: }
2617: } else {
2618: for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
2619: }
2620: /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
2621: for (i = 0; i < dimR; i++) {
2622: PetscReal *swap;
2624: for (j = 0; j < (numV / 2); j++) {
2625: for (k = 0; k < dimC; k++) {
2626: cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
2627: cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
2628: }
2629: }
2631: if (i < dimR - 1) {
2632: swap = cellCoeffs;
2633: cellCoeffs = cellCoords;
2634: cellCoords = swap;
2635: }
2636: }
2637: PetscMemzero(refCoords,numPoints * dimR * sizeof (PetscReal));
2638: for (j = 0; j < numPoints; j++) {
2639: for (i = 0; i < maxIts; i++) {
2640: PetscReal *guess = &refCoords[dimR * j];
2642: /* compute -residual and Jacobian */
2643: for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];}
2644: for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;}
2645: for (k = 0; k < numV; k++) {
2646: PetscReal extCoord = 1.;
2647: for (l = 0; l < dimR; l++) {
2648: PetscReal coord = guess[l];
2649: PetscInt dep = (k & (1 << l)) >> l;
2651: extCoord *= dep * coord + !dep;
2652: extJ[l] = dep;
2654: for (m = 0; m < dimR; m++) {
2655: PetscReal coord = guess[m];
2656: PetscInt dep = ((k & (1 << m)) >> m) && (m != l);
2657: PetscReal mult = dep * coord + !dep;
2659: extJ[l] *= mult;
2660: }
2661: }
2662: for (l = 0; l < dimC; l++) {
2663: PetscReal coeff = cellCoeffs[dimC * k + l];
2665: resNeg[l] -= coeff * extCoord;
2666: for (m = 0; m < dimR; m++) {
2667: J[dimR * l + m] += coeff * extJ[m];
2668: }
2669: }
2670: }
2671: #if 0 && defined(PETSC_USE_DEBUG)
2672: {
2673: PetscReal maxAbs = 0.;
2675: for (l = 0; l < dimC; l++) {
2676: maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
2677: }
2678: PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);
2679: }
2680: #endif
2682: DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);
2683: }
2684: }
2685: DMRestoreWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);
2686: DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);
2687: DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);
2688: return(0);
2689: }
2691: static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR)
2692: {
2693: PetscInt coordSize, i, j, k, l, numV = (1 << dimR);
2694: PetscScalar *coordsScalar = NULL;
2695: PetscReal *cellData, *cellCoords, *cellCoeffs;
2700: DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);
2701: if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);
2702: DMGetWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);
2703: cellCoords = &cellData[0];
2704: cellCoeffs = &cellData[coordSize];
2705: if (dimR == 2) {
2706: const PetscInt zToPlex[4] = {0, 1, 3, 2};
2708: for (i = 0; i < 4; i++) {
2709: PetscInt plexI = zToPlex[i];
2711: for (j = 0; j < dimC; j++) {
2712: cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2713: }
2714: }
2715: } else if (dimR == 3) {
2716: const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6};
2718: for (i = 0; i < 8; i++) {
2719: PetscInt plexI = zToPlex[i];
2721: for (j = 0; j < dimC; j++) {
2722: cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]);
2723: }
2724: }
2725: } else {
2726: for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);}
2727: }
2728: /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */
2729: for (i = 0; i < dimR; i++) {
2730: PetscReal *swap;
2732: for (j = 0; j < (numV / 2); j++) {
2733: for (k = 0; k < dimC; k++) {
2734: cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]);
2735: cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]);
2736: }
2737: }
2739: if (i < dimR - 1) {
2740: swap = cellCoeffs;
2741: cellCoeffs = cellCoords;
2742: cellCoords = swap;
2743: }
2744: }
2745: PetscMemzero(realCoords,numPoints * dimC * sizeof (PetscReal));
2746: for (j = 0; j < numPoints; j++) {
2747: const PetscReal *guess = &refCoords[dimR * j];
2748: PetscReal *mapped = &realCoords[dimC * j];
2750: for (k = 0; k < numV; k++) {
2751: PetscReal extCoord = 1.;
2752: for (l = 0; l < dimR; l++) {
2753: PetscReal coord = guess[l];
2754: PetscInt dep = (k & (1 << l)) >> l;
2756: extCoord *= dep * coord + !dep;
2757: }
2758: for (l = 0; l < dimC; l++) {
2759: PetscReal coeff = cellCoeffs[dimC * k + l];
2761: mapped[l] += coeff * extCoord;
2762: }
2763: }
2764: }
2765: DMRestoreWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);
2766: DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);
2767: return(0);
2768: }
2770: /* TODO: TOBY please fix this for Nc > 1 */
2771: static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
2772: {
2773: PetscInt numComp, pdim, i, j, k, l, m, maxIter = 7, coordSize;
2774: PetscScalar *nodes = NULL;
2775: PetscReal *invV, *modes;
2776: PetscReal *B, *D, *resNeg;
2777: PetscScalar *J, *invJ, *work;
2781: PetscFEGetDimension(fe, &pdim);
2782: PetscFEGetNumComponents(fe, &numComp);
2783: if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
2784: DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);
2785: /* convert nodes to values in the stable evaluation basis */
2786: DMGetWorkArray(dm,pdim,PETSC_REAL,&modes);
2787: invV = fe->invV;
2788: for (i = 0; i < pdim; ++i) {
2789: modes[i] = 0.;
2790: for (j = 0; j < pdim; ++j) {
2791: modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
2792: }
2793: }
2794: DMGetWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,PETSC_REAL,&B);
2795: D = &B[pdim*Nc];
2796: resNeg = &D[pdim*Nc * dimR];
2797: DMGetWorkArray(dm,3 * Nc * dimR,PETSC_SCALAR,&J);
2798: invJ = &J[Nc * dimR];
2799: work = &invJ[Nc * dimR];
2800: for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;}
2801: for (j = 0; j < numPoints; j++) {
2802: for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */
2803: PetscReal *guess = &refCoords[j * dimR];
2804: PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);
2805: for (k = 0; k < Nc; k++) {resNeg[k] = realCoords[j * Nc + k];}
2806: for (k = 0; k < Nc * dimR; k++) {J[k] = 0.;}
2807: for (k = 0; k < pdim; k++) {
2808: for (l = 0; l < Nc; l++) {
2809: resNeg[l] -= modes[k] * B[k * Nc + l];
2810: for (m = 0; m < dimR; m++) {
2811: J[l * dimR + m] += modes[k] * D[(k * Nc + l) * dimR + m];
2812: }
2813: }
2814: }
2815: #if 0 && defined(PETSC_USE_DEBUG)
2816: {
2817: PetscReal maxAbs = 0.;
2819: for (l = 0; l < Nc; l++) {
2820: maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l]));
2821: }
2822: PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);
2823: }
2824: #endif
2825: DMPlexCoordinatesToReference_NewtonUpdate(Nc,dimR,J,invJ,work,resNeg,guess);
2826: }
2827: }
2828: DMRestoreWorkArray(dm,3 * Nc * dimR,PETSC_SCALAR,&J);
2829: DMRestoreWorkArray(dm,pdim * Nc + pdim * Nc * dimR + Nc,PETSC_REAL,&B);
2830: DMRestoreWorkArray(dm,pdim,PETSC_REAL,&modes);
2831: DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);
2832: return(0);
2833: }
2835: /* TODO: TOBY please fix this for Nc > 1 */
2836: static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt Nc, PetscInt dimR)
2837: {
2838: PetscInt numComp, pdim, i, j, k, l, coordSize;
2839: PetscScalar *nodes = NULL;
2840: PetscReal *invV, *modes;
2841: PetscReal *B;
2845: PetscFEGetDimension(fe, &pdim);
2846: PetscFEGetNumComponents(fe, &numComp);
2847: if (numComp != Nc) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,Nc);
2848: DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);
2849: /* convert nodes to values in the stable evaluation basis */
2850: DMGetWorkArray(dm,pdim,PETSC_REAL,&modes);
2851: invV = fe->invV;
2852: for (i = 0; i < pdim; ++i) {
2853: modes[i] = 0.;
2854: for (j = 0; j < pdim; ++j) {
2855: modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]);
2856: }
2857: }
2858: DMGetWorkArray(dm,numPoints * pdim * Nc,PETSC_REAL,&B);
2859: PetscSpaceEvaluate(fe->basisSpace, numPoints, refCoords, B, NULL, NULL);
2860: for (i = 0; i < numPoints * Nc; i++) {realCoords[i] = 0.;}
2861: for (j = 0; j < numPoints; j++) {
2862: PetscReal *mapped = &realCoords[j * Nc];
2864: for (k = 0; k < pdim; k++) {
2865: for (l = 0; l < Nc; l++) {
2866: mapped[l] += modes[k] * B[(j * pdim + k) * Nc + l];
2867: }
2868: }
2869: }
2870: DMRestoreWorkArray(dm,numPoints * pdim * Nc,PETSC_REAL,&B);
2871: DMRestoreWorkArray(dm,pdim,PETSC_REAL,&modes);
2872: DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);
2873: return(0);
2874: }
2876: /*@
2877: DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element
2878: map. This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not
2879: extend uniquely outside the reference cell (e.g, most non-affine maps)
2881: Not collective
2883: Input Parameters:
2884: + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
2885: implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
2886: as a multilinear map for tensor-product elements
2887: . cell - the cell whose map is used.
2888: . numPoints - the number of points to locate
2889: - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
2891: Output Parameters:
2892: . refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
2894: Level: intermediate
2895: @*/
2896: PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
2897: {
2898: PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
2899: DM coordDM = NULL;
2900: Vec coords;
2901: PetscFE fe = NULL;
2906: DMGetDimension(dm,&dimR);
2907: DMGetCoordinateDim(dm,&dimC);
2908: if (dimR <= 0 || dimC <= 0 || numPoints <= 0) return(0);
2909: DMPlexGetDepth(dm,&depth);
2910: DMGetCoordinatesLocal(dm,&coords);
2911: DMGetCoordinateDM(dm,&coordDM);
2912: if (coordDM) {
2913: PetscInt coordFields;
2915: DMGetNumFields(coordDM,&coordFields);
2916: if (coordFields) {
2917: PetscClassId id;
2918: PetscObject disc;
2920: DMGetField(coordDM,0,&disc);
2921: PetscObjectGetClassId(disc,&id);
2922: if (id == PETSCFE_CLASSID) {
2923: fe = (PetscFE) disc;
2924: }
2925: }
2926: }
2927: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
2928: DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);
2929: cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
2930: if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
2931: if (!fe) { /* implicit discretization: affine or multilinear */
2932: PetscInt coneSize;
2933: PetscBool isSimplex, isTensor;
2935: DMPlexGetConeSize(dm,cell,&coneSize);
2936: isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
2937: isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
2938: if (isSimplex) {
2939: PetscReal detJ, *v0, *J, *invJ;
2941: DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);
2942: J = &v0[dimC];
2943: invJ = &J[dimC * dimC];
2944: DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);
2945: for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
2946: CoordinatesRealToRef(dimC, dimR, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]);
2947: }
2948: DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);
2949: } else if (isTensor) {
2950: DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);
2951: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
2952: } else {
2953: DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);
2954: }
2955: return(0);
2956: }
2958: /*@
2959: DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map.
2961: Not collective
2963: Input Parameters:
2964: + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
2965: implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
2966: as a multilinear map for tensor-product elements
2967: . cell - the cell whose map is used.
2968: . numPoints - the number of points to locate
2969: + refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension())
2971: Output Parameters:
2972: . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())
2974: Level: intermediate
2975: @*/
2976: PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[])
2977: {
2978: PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i;
2979: DM coordDM = NULL;
2980: Vec coords;
2981: PetscFE fe = NULL;
2986: DMGetDimension(dm,&dimR);
2987: DMGetCoordinateDim(dm,&dimC);
2988: if (dimR <= 0 || dimC <= 0 || numPoints <= 0) return(0);
2989: DMPlexGetDepth(dm,&depth);
2990: DMGetCoordinatesLocal(dm,&coords);
2991: DMGetCoordinateDM(dm,&coordDM);
2992: if (coordDM) {
2993: PetscInt coordFields;
2995: DMGetNumFields(coordDM,&coordFields);
2996: if (coordFields) {
2997: PetscClassId id;
2998: PetscObject disc;
3000: DMGetField(coordDM,0,&disc);
3001: PetscObjectGetClassId(disc,&id);
3002: if (id == PETSCFE_CLASSID) {
3003: fe = (PetscFE) disc;
3004: }
3005: }
3006: }
3007: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
3008: DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);
3009: cEnd = cEndInterior > 0 ? cEndInterior : cEnd;
3010: if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);
3011: if (!fe) { /* implicit discretization: affine or multilinear */
3012: PetscInt coneSize;
3013: PetscBool isSimplex, isTensor;
3015: DMPlexGetConeSize(dm,cell,&coneSize);
3016: isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE;
3017: isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE;
3018: if (isSimplex) {
3019: PetscReal detJ, *v0, *J;
3021: DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);
3022: J = &v0[dimC];
3023: DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);
3024: for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */
3025: CoordinatesRefToReal(dimC, dimR, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]);
3026: }
3027: DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);
3028: } else if (isTensor) {
3029: DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);
3030: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize);
3031: } else {
3032: DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);
3033: }
3034: return(0);
3035: }