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