xref: /petsc/src/dm/impls/plex/tests/ex8.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Tests for cell geometry\n\n";
2 
3 #include <petscdmplex.h>
4 
5 typedef enum {RUN_REFERENCE, RUN_HEX_CURVED, RUN_FILE, RUN_DISPLAY} RunType;
6 
7 typedef struct {
8   DM        dm;
9   RunType   runType;                      /* Type of mesh to use */
10   PetscBool transform;                    /* Use random coordinate transformations */
11   /* Data for input meshes */
12   PetscReal *v0, *J, *invJ, *detJ;        /* FEM data */
13   PetscReal *centroid, *normal, *vol;     /* FVM data */
14 } AppCtx;
15 
16 static PetscErrorCode ReadMesh(MPI_Comm comm, AppCtx *user, DM *dm)
17 {
18   PetscFunctionBegin;
19   PetscCall(DMCreate(comm, dm));
20   PetscCall(DMSetType(*dm, DMPLEX));
21   PetscCall(DMSetFromOptions(*dm));
22   PetscCall(DMSetApplicationContext(*dm, user));
23   PetscCall(PetscObjectSetName((PetscObject) *dm, "Input Mesh"));
24   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
25   PetscFunctionReturn(0);
26 }
27 
28 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
29 {
30   const char    *runTypes[4] = {"reference", "hex_curved", "file", "display"};
31   PetscInt       run;
32 
33   PetscFunctionBeginUser;
34   options->runType   = RUN_REFERENCE;
35   options->transform = PETSC_FALSE;
36 
37   PetscOptionsBegin(comm, "", "Geometry Test Options", "DMPLEX");
38   run  = options->runType;
39   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex8.c", runTypes, 3, runTypes[options->runType], &run, NULL));
40   options->runType = (RunType) run;
41   PetscCall(PetscOptionsBool("-transform", "Use random transforms", "ex8.c", options->transform, &options->transform, NULL));
42 
43   if (options->runType == RUN_FILE) {
44     PetscInt  dim, cStart, cEnd, numCells, n;
45     PetscBool flg, feFlg;
46 
47     PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm));
48     PetscCall(DMGetDimension(options->dm, &dim));
49     PetscCall(DMPlexGetHeightStratum(options->dm, 0, &cStart, &cEnd));
50     numCells = cEnd-cStart;
51     PetscCall(PetscMalloc4(numCells*dim, &options->v0, numCells*dim*dim, &options->J, numCells*dim*dim, &options->invJ, numCells, &options->detJ));
52     PetscCall(PetscMalloc1(numCells*dim, &options->centroid));
53     PetscCall(PetscMalloc1(numCells*dim, &options->normal));
54     PetscCall(PetscMalloc1(numCells, &options->vol));
55     n = numCells*dim;
56     PetscCall(PetscOptionsRealArray("-v0", "Input v0 for each cell", "ex8.c", options->v0, &n, &feFlg));
57     PetscCheck(!feFlg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of v0 %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim);
58     n = numCells*dim*dim;
59     PetscCall(PetscOptionsRealArray("-J", "Input Jacobian for each cell", "ex8.c", options->J, &n, &flg));
60     PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of J %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim*dim);
61     n = numCells*dim*dim;
62     PetscCall(PetscOptionsRealArray("-invJ", "Input inverse Jacobian for each cell", "ex8.c", options->invJ, &n, &flg));
63     PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of invJ %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim*dim);
64     n = numCells;
65     PetscCall(PetscOptionsRealArray("-detJ", "Input Jacobian determinant for each cell", "ex8.c", options->detJ, &n, &flg));
66     PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of detJ %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells);
67     n = numCells*dim;
68     if (!feFlg) {
69       PetscCall(PetscFree4(options->v0, options->J, options->invJ, options->detJ));
70       options->v0 = options->J = options->invJ = options->detJ = NULL;
71     }
72     PetscCall(PetscOptionsRealArray("-centroid", "Input centroid for each cell", "ex8.c", options->centroid, &n, &flg));
73     PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of centroid %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim);
74     if (!flg) {
75       PetscCall(PetscFree(options->centroid));
76       options->centroid = NULL;
77     }
78     n = numCells*dim;
79     PetscCall(PetscOptionsRealArray("-normal", "Input normal for each cell", "ex8.c", options->normal, &n, &flg));
80     PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of normal %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim);
81     if (!flg) {
82       PetscCall(PetscFree(options->normal));
83       options->normal = NULL;
84     }
85     n = numCells;
86     PetscCall(PetscOptionsRealArray("-vol", "Input volume for each cell", "ex8.c", options->vol, &n, &flg));
87     PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of vol %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells);
88     if (!flg) {
89       PetscCall(PetscFree(options->vol));
90       options->vol = NULL;
91     }
92   } else if (options->runType == RUN_DISPLAY) {
93     PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm));
94   }
95   PetscOptionsEnd();
96 
97   if (options->transform) PetscCall(PetscPrintf(comm, "Using random transforms\n"));
98   PetscFunctionReturn(0);
99 }
100 
101 static PetscErrorCode ChangeCoordinates(DM dm, PetscInt spaceDim, PetscScalar vertexCoords[])
102 {
103   PetscSection   coordSection;
104   Vec            coordinates;
105   PetscScalar   *coords;
106   PetscInt       vStart, vEnd, v, d, coordSize;
107 
108   PetscFunctionBegin;
109   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
110   PetscCall(DMGetCoordinateSection(dm, &coordSection));
111   PetscCall(PetscSectionSetNumFields(coordSection, 1));
112   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim));
113   PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd));
114   for (v = vStart; v < vEnd; ++v) {
115     PetscCall(PetscSectionSetDof(coordSection, v, spaceDim));
116     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim));
117   }
118   PetscCall(PetscSectionSetUp(coordSection));
119   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
120   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
121   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
122   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
123   PetscCall(VecSetFromOptions(coordinates));
124   PetscCall(VecGetArray(coordinates, &coords));
125   for (v = vStart; v < vEnd; ++v) {
126     PetscInt off;
127 
128     PetscCall(PetscSectionGetOffset(coordSection, v, &off));
129     for (d = 0; d < spaceDim; ++d) {
130       coords[off+d] = vertexCoords[(v-vStart)*spaceDim+d];
131     }
132   }
133   PetscCall(VecRestoreArray(coordinates, &coords));
134   PetscCall(DMSetCoordinateDim(dm, spaceDim));
135   PetscCall(DMSetCoordinatesLocal(dm, coordinates));
136   PetscCall(VecDestroy(&coordinates));
137   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
138   PetscFunctionReturn(0);
139 }
140 
141 #define RelativeError(a,b) PetscAbs(a-b)/(1.0+PetscMax(PetscAbs(a),PetscAbs(b)))
142 
143 static PetscErrorCode CheckFEMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx)
144 {
145   PetscReal      v0[3], J[9], invJ[9], detJ;
146   PetscInt       d, i, j;
147 
148   PetscFunctionBegin;
149   PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ));
150   for (d = 0; d < spaceDim; ++d) {
151     if (v0[d] != v0Ex[d]) {
152       switch (spaceDim) {
153       case 2: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid v0 (%g, %g) != (%g, %g)", (double)v0[0], (double)v0[1], (double)v0Ex[0], (double)v0Ex[1]);
154       case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid v0 (%g, %g, %g) != (%g, %g, %g)", (double)v0[0], (double)v0[1], (double)v0[2], (double)v0Ex[0], (double)v0Ex[1], (double)v0Ex[2]);
155       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid space dimension %" PetscInt_FMT, spaceDim);
156       }
157     }
158   }
159   for (i = 0; i < spaceDim; ++i) {
160     for (j = 0; j < spaceDim; ++j) {
161       PetscCheck(RelativeError(J[i*spaceDim+j], JEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid J[%" PetscInt_FMT ",%" PetscInt_FMT "]: %g != %g", i, j, (double)J[i*spaceDim+j], (double)JEx[i*spaceDim+j]);
162       PetscCheck(RelativeError(invJ[i*spaceDim+j], invJEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid invJ[%" PetscInt_FMT ",%" PetscInt_FMT "]: %g != %g", i, j, (double)invJ[i*spaceDim+j], (double)invJEx[i*spaceDim+j]);
163     }
164   }
165   PetscCheck(RelativeError(detJ, detJEx) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid |J| = %g != %g diff %g", (double)detJ, (double)detJEx,(double)(detJ - detJEx));
166   PetscFunctionReturn(0);
167 }
168 
169 static PetscErrorCode CheckFVMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx)
170 {
171   PetscReal      tol = PetscMax(10*PETSC_SMALL, 1e-10);
172   PetscReal      centroid[3], normal[3], vol;
173   PetscInt       d;
174 
175   PetscFunctionBegin;
176   PetscCall(DMPlexComputeCellGeometryFVM(dm, cell, volEx? &vol : NULL, centroidEx? centroid : NULL, normalEx? normal : NULL));
177   for (d = 0; d < spaceDim; ++d) {
178     if (centroidEx)
179       PetscCheck(RelativeError(centroid[d], centroidEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid centroid[%" PetscInt_FMT "]: %g != %g diff %g", cell, d, (double)centroid[d], (double)centroidEx[d],(double)(centroid[d]-centroidEx[d]));
180     if (normalEx)
181       PetscCheck(RelativeError(normal[d], normalEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid normal[%" PetscInt_FMT "]: %g != %g", cell, d, (double) normal[d], (double) normalEx[d]);
182   }
183   if (volEx)
184     PetscCheck(RelativeError(volEx, vol) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid volume = %g != %g diff %g", cell, (double)vol, (double)volEx,(double)(vol - volEx));
185   PetscFunctionReturn(0);
186 }
187 
188 static PetscErrorCode CheckGaussLaw(DM dm, PetscInt cell)
189 {
190   DMPolytopeType  ct;
191   PetscReal       tol = PetscMax(10*PETSC_SMALL, 1e-10);
192   PetscReal       normal[3], integral[3] = {0., 0., 0.}, area;
193   const PetscInt *cone, *ornt;
194   PetscInt        coneSize, f, dim, cdim, d;
195 
196   PetscFunctionBegin;
197   PetscCall(DMGetDimension(dm, &dim));
198   PetscCall(DMGetCoordinateDim(dm, &cdim));
199   if (dim != cdim) PetscFunctionReturn(0);
200   PetscCall(DMPlexGetCellType(dm, cell, &ct));
201   if (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) PetscFunctionReturn(0);
202   PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
203   PetscCall(DMPlexGetCone(dm, cell, &cone));
204   PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt));
205   for (f = 0; f < coneSize; ++f) {
206     const PetscInt sgn = dim == 1? (f == 0 ? -1 : 1) : (ornt[f] < 0 ? -1 : 1);
207 
208     PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], &area, NULL, normal));
209     for (d = 0; d < cdim; ++d) integral[d] += sgn*area*normal[d];
210   }
211   for (d = 0; d < cdim; ++d) PetscCheck(PetscAbsReal(integral[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT " Surface integral for component %" PetscInt_FMT ": %g != 0. as it should be for a constant field", cell, d, (double) integral[d]);
212   PetscFunctionReturn(0);
213 }
214 
215 static PetscErrorCode CheckCell(DM dm, PetscInt cell, PetscBool transform, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx, PetscReal faceCentroidEx[], PetscReal faceNormalEx[], PetscReal faceVolEx[])
216 {
217   const PetscInt *cone;
218   PetscInt        coneSize, c;
219   PetscInt        dim, depth, cdim;
220 
221   PetscFunctionBegin;
222   PetscCall(DMPlexGetDepth(dm, &depth));
223   PetscCall(DMGetDimension(dm, &dim));
224   PetscCall(DMGetCoordinateDim(dm, &cdim));
225   if (v0Ex) PetscCall(CheckFEMGeometry(dm, cell, cdim, v0Ex, JEx, invJEx, detJEx));
226   if (dim == depth && centroidEx) {
227     PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidEx, normalEx, volEx));
228     PetscCall(CheckGaussLaw(dm, cell));
229     if (faceCentroidEx) {
230       PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
231       PetscCall(DMPlexGetCone(dm, cell, &cone));
232       for (c = 0; c < coneSize; ++c) {
233         PetscCall(CheckFVMGeometry(dm, cone[c], dim, &faceCentroidEx[c*dim], &faceNormalEx[c*dim], faceVolEx[c]));
234       }
235     }
236   }
237   if (transform) {
238     Vec          coordinates;
239     PetscSection coordSection;
240     PetscScalar *coords = NULL, *origCoords, *newCoords;
241     PetscReal   *v0ExT, *JExT, *invJExT, detJExT=0, *centroidExT, *normalExT, volExT=0;
242     PetscReal   *faceCentroidExT, *faceNormalExT, faceVolExT;
243     PetscRandom  r, ang, ang2;
244     PetscInt     coordSize, numCorners, t;
245 
246     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
247     PetscCall(DMGetCoordinateSection(dm, &coordSection));
248     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords));
249     PetscCall(PetscMalloc2(coordSize, &origCoords, coordSize, &newCoords));
250     PetscCall(PetscMalloc5(cdim, &v0ExT, cdim*cdim, &JExT, cdim*cdim, &invJExT, cdim, &centroidExT, cdim, &normalExT));
251     PetscCall(PetscMalloc2(cdim, &faceCentroidExT, cdim, &faceNormalExT));
252     for (c = 0; c < coordSize; ++c) origCoords[c] = coords[c];
253     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords));
254     numCorners = coordSize/cdim;
255 
256     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
257     PetscCall(PetscRandomSetFromOptions(r));
258     PetscCall(PetscRandomSetInterval(r, 0.0, 10.0));
259     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang));
260     PetscCall(PetscRandomSetFromOptions(ang));
261     PetscCall(PetscRandomSetInterval(ang, 0.0, 2*PETSC_PI));
262     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang2));
263     PetscCall(PetscRandomSetFromOptions(ang2));
264     PetscCall(PetscRandomSetInterval(ang2, 0.0, PETSC_PI));
265     for (t = 0; t < 1; ++t) {
266       PetscScalar trans[3];
267       PetscReal   R[9], rot[3], rotM[9];
268       PetscReal   scale, phi, theta, psi = 0.0, norm;
269       PetscInt    d, e, f, p;
270 
271       for (c = 0; c < coordSize; ++c) newCoords[c] = origCoords[c];
272       PetscCall(PetscRandomGetValueReal(r, &scale));
273       PetscCall(PetscRandomGetValueReal(ang, &phi));
274       PetscCall(PetscRandomGetValueReal(ang2, &theta));
275       for (d = 0; d < cdim; ++d) {
276         PetscCall(PetscRandomGetValue(r, &trans[d]));
277       }
278       switch (cdim) {
279       case 2:
280         R[0] = PetscCosReal(phi); R[1] = -PetscSinReal(phi);
281         R[2] = PetscSinReal(phi); R[3] =  PetscCosReal(phi);
282         break;
283       case 3:
284       {
285         const PetscReal ct = PetscCosReal(theta), st = PetscSinReal(theta);
286         const PetscReal cp = PetscCosReal(phi),   sp = PetscSinReal(phi);
287         const PetscReal cs = PetscCosReal(psi),   ss = PetscSinReal(psi);
288         R[0] = ct*cs; R[1] = sp*st*cs - cp*ss;    R[2] = sp*ss    + cp*st*cs;
289         R[3] = ct*ss; R[4] = cp*cs    + sp*st*ss; R[5] = cp*st*ss - sp*cs;
290         R[6] = -st;   R[7] = sp*ct;               R[8] = cp*ct;
291         break;
292       }
293       default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid coordinate dimension %" PetscInt_FMT, cdim);
294       }
295       if (v0Ex) {
296         detJExT = detJEx;
297         for (d = 0; d < cdim; ++d) {
298           v0ExT[d] = v0Ex[d];
299           for (e = 0; e < cdim; ++e) {
300             JExT[d*cdim+e]    = JEx[d*cdim+e];
301             invJExT[d*cdim+e] = invJEx[d*cdim+e];
302           }
303         }
304         for (d = 0; d < cdim; ++d) {
305           v0ExT[d] *= scale;
306           v0ExT[d] += PetscRealPart(trans[d]);
307           /* Only scale dimensions in the manifold */
308           for (e = 0; e < dim; ++e) {
309             JExT[d*cdim+e]    *= scale;
310             invJExT[d*cdim+e] /= scale;
311           }
312           if (d < dim) detJExT *= scale;
313         }
314         /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */
315         for (d = 0; d < cdim; ++d) {
316           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
317             rot[d] += R[d*cdim+e] * v0ExT[e];
318           }
319         }
320         for (d = 0; d < cdim; ++d) v0ExT[d] = rot[d];
321         for (d = 0; d < cdim; ++d) {
322           for (e = 0; e < cdim; ++e) {
323             for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) {
324               rotM[d*cdim+e] += R[d*cdim+f] * JExT[f*cdim+e];
325             }
326           }
327         }
328         for (d = 0; d < cdim; ++d) {
329           for (e = 0; e < cdim; ++e) {
330             JExT[d*cdim+e] = rotM[d*cdim+e];
331           }
332         }
333         for (d = 0; d < cdim; ++d) {
334           for (e = 0; e < cdim; ++e) {
335             for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) {
336               rotM[d*cdim+e] += invJExT[d*cdim+f] * R[e*cdim+f];
337             }
338           }
339         }
340         for (d = 0; d < cdim; ++d) {
341           for (e = 0; e < cdim; ++e) {
342             invJExT[d*cdim+e] = rotM[d*cdim+e];
343           }
344         }
345       }
346       if (centroidEx) {
347         volExT = volEx;
348         for (d = 0; d < cdim; ++d) {
349           centroidExT[d]  = centroidEx[d];
350           normalExT[d]    = normalEx[d];
351         }
352         for (d = 0; d < cdim; ++d) {
353           centroidExT[d] *= scale;
354           centroidExT[d] += PetscRealPart(trans[d]);
355           normalExT[d]   /= scale;
356           /* Only scale dimensions in the manifold */
357           if (d < dim) volExT  *= scale;
358         }
359         /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */
360         for (d = 0; d < cdim; ++d) {
361           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
362             rot[d] += R[d*cdim+e] * centroidExT[e];
363           }
364         }
365         for (d = 0; d < cdim; ++d) centroidExT[d] = rot[d];
366         for (d = 0; d < cdim; ++d) {
367           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
368             rot[d] += R[d*cdim+e] * normalExT[e];
369           }
370         }
371         for (d = 0; d < cdim; ++d) normalExT[d] = rot[d];
372         for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(normalExT[d]);
373         norm = PetscSqrtReal(norm);
374         if (norm != 0.) for (d = 0; d < cdim; ++d) normalExT[d] /= norm;
375       }
376       for (d = 0; d < cdim; ++d) {
377         for (p = 0; p < numCorners; ++p) {
378           newCoords[p*cdim+d] *= scale;
379           newCoords[p*cdim+d] += trans[d];
380         }
381       }
382       for (p = 0; p < numCorners; ++p) {
383         for (d = 0; d < cdim; ++d) {
384           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
385             rot[d] += R[d*cdim+e] * PetscRealPart(newCoords[p*cdim+e]);
386           }
387         }
388         for (d = 0; d < cdim; ++d) newCoords[p*cdim+d] = rot[d];
389       }
390 
391       PetscCall(ChangeCoordinates(dm, cdim, newCoords));
392       if (v0Ex) PetscCall(CheckFEMGeometry(dm, 0, cdim, v0ExT, JExT, invJExT, detJExT));
393       if (dim == depth && centroidEx) {
394         PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidExT, normalExT, volExT));
395         PetscCall(CheckGaussLaw(dm, cell));
396         if (faceCentroidEx) {
397           PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
398           PetscCall(DMPlexGetCone(dm, cell, &cone));
399           for (c = 0; c < coneSize; ++c) {
400             PetscInt off = c*cdim;
401 
402             faceVolExT = faceVolEx[c];
403             for (d = 0; d < cdim; ++d) {
404               faceCentroidExT[d]  = faceCentroidEx[off+d];
405               faceNormalExT[d]    = faceNormalEx[off+d];
406             }
407             for (d = 0; d < cdim; ++d) {
408               faceCentroidExT[d] *= scale;
409               faceCentroidExT[d] += PetscRealPart(trans[d]);
410               faceNormalExT[d]   /= scale;
411               /* Only scale dimensions in the manifold */
412               if (d < dim-1) {
413                 faceVolExT *= scale;
414               }
415             }
416             for (d = 0; d < cdim; ++d) {
417               for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
418                 rot[d] += R[d*cdim+e] * faceCentroidExT[e];
419               }
420             }
421             for (d = 0; d < cdim; ++d) faceCentroidExT[d] = rot[d];
422             for (d = 0; d < cdim; ++d) {
423               for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
424                 rot[d] += R[d*cdim+e] * faceNormalExT[e];
425               }
426             }
427             for (d = 0; d < cdim; ++d) faceNormalExT[d] = rot[d];
428             for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(faceNormalExT[d]);
429             norm = PetscSqrtReal(norm);
430             for (d = 0; d < cdim; ++d) faceNormalExT[d] /= norm;
431 
432             PetscCall(CheckFVMGeometry(dm, cone[c], cdim, faceCentroidExT, faceNormalExT, faceVolExT));
433           }
434         }
435       }
436     }
437     PetscCall(PetscRandomDestroy(&r));
438     PetscCall(PetscRandomDestroy(&ang));
439     PetscCall(PetscRandomDestroy(&ang2));
440     PetscCall(PetscFree2(origCoords, newCoords));
441     PetscCall(PetscFree5(v0ExT, JExT, invJExT, centroidExT, normalExT));
442     PetscCall(PetscFree2(faceCentroidExT, faceNormalExT));
443   }
444   PetscFunctionReturn(0);
445 }
446 
447 static PetscErrorCode TestTriangle(MPI_Comm comm, PetscBool transform)
448 {
449   DM             dm;
450 
451   PetscFunctionBegin;
452   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRIANGLE, &dm));
453   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
454   /* Check reference geometry: determinant is scaled by reference volume (2.0) */
455   {
456     PetscReal v0Ex[2]       = {-1.0, -1.0};
457     PetscReal JEx[4]        = {1.0, 0.0, 0.0, 1.0};
458     PetscReal invJEx[4]     = {1.0, 0.0, 0.0, 1.0};
459     PetscReal detJEx        = 1.0;
460     PetscReal centroidEx[2] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.)};
461     PetscReal normalEx[2]   = {0.0, 0.0};
462     PetscReal volEx         = 2.0;
463 
464     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
465   }
466   /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (2.0) */
467   {
468     PetscScalar vertexCoords[9] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, -1.0, 1.0, 0.0};
469     PetscReal   v0Ex[3]         = {-1.0, -1.0, 0.0};
470     PetscReal   JEx[9]          = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
471     PetscReal   invJEx[9]       = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
472     PetscReal   detJEx          = 1.0;
473     PetscReal   centroidEx[3]   = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0};
474     PetscReal   normalEx[3]     = {0.0, 0.0, 1.0};
475     PetscReal   volEx           = 2.0;
476 
477     PetscCall(ChangeCoordinates(dm, 3, vertexCoords));
478     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
479   }
480   /* Cleanup */
481   PetscCall(DMDestroy(&dm));
482   PetscFunctionReturn(0);
483 }
484 
485 static PetscErrorCode TestQuadrilateral(MPI_Comm comm, PetscBool transform)
486 {
487   DM             dm;
488 
489   PetscFunctionBegin;
490   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_QUADRILATERAL, &dm));
491   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
492   /* Check reference geometry: determinant is scaled by reference volume (2.0) */
493   {
494     PetscReal v0Ex[2]       = {-1.0, -1.0};
495     PetscReal JEx[4]        = {1.0, 0.0, 0.0, 1.0};
496     PetscReal invJEx[4]     = {1.0, 0.0, 0.0, 1.0};
497     PetscReal detJEx        = 1.0;
498     PetscReal centroidEx[2] = {0.0, 0.0};
499     PetscReal normalEx[2]   = {0.0, 0.0};
500     PetscReal volEx         = 4.0;
501 
502     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
503   }
504   /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (4.0) */
505   {
506     PetscScalar vertexCoords[12] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, 1.0, 1.0, 0.0, -1.0, 1.0, 0.0};
507     PetscReal   v0Ex[3]          = {-1.0, -1.0, 0.0};
508     PetscReal   JEx[9]           = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
509     PetscReal   invJEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
510     PetscReal   detJEx           = 1.0;
511     PetscReal   centroidEx[3]    = {0.0, 0.0, 0.0};
512     PetscReal   normalEx[3]      = {0.0, 0.0, 1.0};
513     PetscReal   volEx            = 4.0;
514 
515     PetscCall(ChangeCoordinates(dm, 3, vertexCoords));
516     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
517   }
518   /* Cleanup */
519   PetscCall(DMDestroy(&dm));
520   PetscFunctionReturn(0);
521 }
522 
523 static PetscErrorCode TestTetrahedron(MPI_Comm comm, PetscBool transform)
524 {
525   DM             dm;
526 
527   PetscFunctionBegin;
528   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TETRAHEDRON, &dm));
529   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
530   /* Check reference geometry: determinant is scaled by reference volume (4/3) */
531   {
532     PetscReal v0Ex[3]       = {-1.0, -1.0, -1.0};
533     PetscReal JEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
534     PetscReal invJEx[9]     = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
535     PetscReal detJEx        = 1.0;
536     PetscReal centroidEx[3] = {-0.5, -0.5, -0.5};
537     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
538     PetscReal volEx         = (PetscReal)4.0/(PetscReal)3.0;
539 
540     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
541   }
542   /* Cleanup */
543   PetscCall(DMDestroy(&dm));
544   PetscFunctionReturn(0);
545 }
546 
547 static PetscErrorCode TestHexahedron(MPI_Comm comm, PetscBool transform)
548 {
549   DM             dm;
550 
551   PetscFunctionBegin;
552   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm));
553   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
554   /* Check reference geometry: determinant is scaled by reference volume 8.0 */
555   {
556     PetscReal v0Ex[3]       = {-1.0, -1.0, -1.0};
557     PetscReal JEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
558     PetscReal invJEx[9]     = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
559     PetscReal detJEx        = 1.0;
560     PetscReal centroidEx[3] = {0.0, 0.0, 0.0};
561     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
562     PetscReal volEx         = 8.0;
563 
564     PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL));
565   }
566   /* Cleanup */
567   PetscCall(DMDestroy(&dm));
568   PetscFunctionReturn(0);
569 }
570 
571 static PetscErrorCode TestHexahedronCurved(MPI_Comm comm)
572 {
573   DM             dm;
574   PetscScalar    coords[24] = {-1.0, -1.0, -1.0,  -1.0,  1.0, -1.0,  1.0, 1.0, -1.0,   1.0, -1.0, -1.0,
575                                -1.0, -1.0,  1.1,   1.0, -1.0,  1.0,  1.0, 1.0,  1.1,  -1.0,  1.0,  1.0};
576 
577   PetscFunctionBegin;
578   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm));
579   PetscCall(ChangeCoordinates(dm, 3, coords));
580   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
581   {
582     PetscReal centroidEx[3] = {0.0, 0.0, 0.016803278688524603};
583     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
584     PetscReal volEx         = 8.1333333333333346;
585 
586     PetscCall(CheckCell(dm, 0, PETSC_FALSE, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, NULL, NULL, NULL));
587   }
588   PetscCall(DMDestroy(&dm));
589   PetscFunctionReturn(0);
590 }
591 
592 /* This wedge is a tensor product cell, rather than a normal wedge */
593 static PetscErrorCode TestWedge(MPI_Comm comm, PetscBool transform)
594 {
595   DM             dm;
596 
597   PetscFunctionBegin;
598   PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRI_PRISM_TENSOR, &dm));
599   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
600   /* Check reference geometry: determinant is scaled by reference volume 4.0 */
601   {
602 #if 0
603     /* FEM geometry is not functional for wedges */
604     PetscReal v0Ex[3]   = {-1.0, -1.0, -1.0};
605     PetscReal JEx[9]    = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
606     PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
607     PetscReal detJEx    = 1.0;
608 #endif
609 
610     {
611       PetscReal       centroidEx[3]      = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0};
612       PetscReal       normalEx[3]        = {0.0, 0.0, 0.0};
613       PetscReal       volEx              = 4.0;
614       PetscReal       faceVolEx[5]       = {2.0, 2.0, 4.0, PETSC_SQRT2*4.0, 4.0};
615       PetscReal       faceNormalEx[15]   = {0.0, 0.0, 1.0,  0.0, 0.0, 1.0,  0.0, -1.0, 0.0,  PETSC_SQRT2/2.0, PETSC_SQRT2/2.0, 0.0,  -1.0, 0.0, 0.0};
616       PetscReal       faceCentroidEx[15] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), -1.0,
617                                             -((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.),  1.0,
618                                             0.0, -1.0, 0.0,  0.0, 0.0, 0.0,  -1.0, 0.0, 0.0};
619 
620       PetscCall(CheckCell(dm, 0, transform, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, faceCentroidEx, faceNormalEx, faceVolEx));
621     }
622   }
623   /* Cleanup */
624   PetscCall(DMDestroy(&dm));
625   PetscFunctionReturn(0);
626 }
627 
628 int main(int argc, char **argv)
629 {
630   AppCtx         user;
631 
632   PetscFunctionBeginUser;
633   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
634   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
635   if (user.runType == RUN_REFERENCE) {
636     PetscCall(TestTriangle(PETSC_COMM_SELF, user.transform));
637     PetscCall(TestQuadrilateral(PETSC_COMM_SELF, user.transform));
638     PetscCall(TestTetrahedron(PETSC_COMM_SELF, user.transform));
639     PetscCall(TestHexahedron(PETSC_COMM_SELF, user.transform));
640     PetscCall(TestWedge(PETSC_COMM_SELF, user.transform));
641   } else if (user.runType == RUN_HEX_CURVED) {
642     PetscCall(TestHexahedronCurved(PETSC_COMM_SELF));
643   } else if (user.runType == RUN_FILE) {
644     PetscInt dim, cStart, cEnd, c;
645 
646     PetscCall(DMGetDimension(user.dm, &dim));
647     PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd));
648     for (c = 0; c < cEnd-cStart; ++c) {
649       PetscReal *v0       = user.v0       ? &user.v0[c*dim] : NULL;
650       PetscReal *J        = user.J        ? &user.J[c*dim*dim] : NULL;
651       PetscReal *invJ     = user.invJ     ? &user.invJ[c*dim*dim] : NULL;
652       PetscReal  detJ     = user.detJ     ?  user.detJ[c] : 0.0;
653       PetscReal *centroid = user.centroid ? &user.centroid[c*dim] : NULL;
654       PetscReal *normal   = user.normal   ? &user.normal[c*dim] : NULL;
655       PetscReal  vol      = user.vol      ?  user.vol[c] : 0.0;
656 
657       PetscCall(CheckCell(user.dm, c+cStart, PETSC_FALSE, v0, J, invJ, detJ, centroid, normal, vol, NULL, NULL, NULL));
658     }
659     PetscCall(PetscFree4(user.v0,user.J,user.invJ,user.detJ));
660     PetscCall(PetscFree(user.centroid));
661     PetscCall(PetscFree(user.normal));
662     PetscCall(PetscFree(user.vol));
663     PetscCall(DMDestroy(&user.dm));
664   } else if (user.runType == RUN_DISPLAY) {
665     DM                 gdm, dmCell;
666     Vec                cellgeom, facegeom;
667     const PetscScalar *cgeom;
668     PetscInt           dim, d, cStart, cEnd, cEndInterior, c;
669 
670     PetscCall(DMGetCoordinateDim(user.dm, &dim));
671     PetscCall(DMPlexConstructGhostCells(user.dm, NULL, NULL, &gdm));
672     if (gdm) {
673       PetscCall(DMDestroy(&user.dm));
674       user.dm = gdm;
675     }
676     PetscCall(DMPlexComputeGeometryFVM(user.dm, &cellgeom, &facegeom));
677     PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd));
678     PetscCall(DMPlexGetGhostCellStratum(user.dm, &cEndInterior, NULL));
679     if (cEndInterior >= 0) cEnd = cEndInterior;
680     PetscCall(VecGetDM(cellgeom, &dmCell));
681     PetscCall(VecGetArrayRead(cellgeom, &cgeom));
682     for (c = 0; c < cEnd-cStart; ++c) {
683       PetscFVCellGeom *cg;
684 
685       PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
686       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %4" PetscInt_FMT ": Centroid (", c));
687       for (d = 0; d < dim; ++d) {
688         if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
689         PetscCall(PetscPrintf(PETSC_COMM_SELF, "%12.2g", (double)cg->centroid[d]));
690       }
691       PetscCall(PetscPrintf(PETSC_COMM_SELF, ") Vol %12.2g\n", (double)cg->volume));
692     }
693     PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
694     PetscCall(VecDestroy(&cellgeom));
695     PetscCall(VecDestroy(&facegeom));
696     PetscCall(DMDestroy(&user.dm));
697   }
698   PetscCall(PetscFinalize());
699   return 0;
700 }
701 
702 /*TEST
703 
704   test:
705     suffix: 1
706     args: -dm_view ascii::ascii_info_detail
707   test:
708     suffix: 2
709     args: -run_type hex_curved
710   test:
711     suffix: 3
712     args: -transform
713   test:
714     suffix: 4
715     requires: exodusii
716     args: -run_type file -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/simpleblock-100.exo -dm_view ascii::ascii_info_detail -v0 -1.5,-0.5,0.5,-0.5,-0.5,0.5,0.5,-0.5,0.5 -J 0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0,0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0,0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0 -invJ 0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0,0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0,0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0 -detJ 0.125,0.125,0.125 -centroid -1.0,0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0 -normal 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0 -vol 1.0,1.0,1.0
717   test:
718     suffix: 5
719     args: -run_type file -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 3,1,1 -dm_plex_box_lower -1.5,-0.5,-0.5 -dm_plex_box_upper 1.5,0.5,0.5 -dm_view ascii::ascii_info_detail -centroid -1.0,0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0 -normal 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0 -vol 1.0,1.0,1.0
720   test:
721     suffix: 6
722     args: -run_type file -dm_plex_dim 1 -dm_plex_simplex 0 -dm_plex_box_faces 3 -dm_plex_box_lower -1.5 -dm_plex_box_upper 1.5 -dm_view ascii::ascii_info_detail -centroid -1.0,0.0,1.0 -vol 1.0,1.0,1.0
723 TEST*/
724