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