xref: /petsc/src/dm/impls/plex/tests/ex8.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
1c4762a1bSJed Brown static char help[] = "Tests for cell geometry\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
54f99dae5SMatthew 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 
164f99dae5SMatthew G. Knepley static PetscErrorCode ReadMesh(MPI_Comm comm, AppCtx *user, DM *dm)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   PetscErrorCode ierr;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBegin;
214f99dae5SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
224f99dae5SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
234f99dae5SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
244f99dae5SMatthew 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 {
324f99dae5SMatthew 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;
489bf2564aSMatt McGurn     PetscBool flg, feFlg;
49c4762a1bSJed Brown 
504f99dae5SMatthew 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;
544f99dae5SMatthew G. Knepley     ierr = PetscMalloc4(numCells*dim, &options->v0, numCells*dim*dim, &options->J, numCells*dim*dim, &options->invJ, numCells, &options->detJ);CHKERRQ(ierr);
559bf2564aSMatt McGurn     ierr = PetscMalloc1(numCells*dim, &options->centroid);CHKERRQ(ierr);
569bf2564aSMatt McGurn     ierr = PetscMalloc1(numCells*dim, &options->normal);CHKERRQ(ierr);
579bf2564aSMatt McGurn     ierr = PetscMalloc1(numCells, &options->vol);CHKERRQ(ierr);
58c4762a1bSJed Brown     n    = numCells*dim;
594f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-v0", "Input v0 for each cell", "ex8.c", options->v0, &n, &feFlg);CHKERRQ(ierr);
60*98921bdaSJacob Faibussowitsch     if (feFlg && n != numCells*dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of v0 %D should be %D", n, numCells*dim);
61c4762a1bSJed Brown     n    = numCells*dim*dim;
624f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-J", "Input Jacobian for each cell", "ex8.c", options->J, &n, &flg);CHKERRQ(ierr);
63*98921bdaSJacob Faibussowitsch     if (flg && n != numCells*dim*dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of J %D should be %D", n, numCells*dim*dim);
64c4762a1bSJed Brown     n    = numCells*dim*dim;
654f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-invJ", "Input inverse Jacobian for each cell", "ex8.c", options->invJ, &n, &flg);CHKERRQ(ierr);
66*98921bdaSJacob Faibussowitsch     if (flg && n != numCells*dim*dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of invJ %D should be %D", n, numCells*dim*dim);
67c4762a1bSJed Brown     n    = numCells;
684f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-detJ", "Input Jacobian determinant for each cell", "ex8.c", options->detJ, &n, &flg);CHKERRQ(ierr);
69*98921bdaSJacob Faibussowitsch     if (flg && n != numCells) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of detJ %D should be %D", n, numCells);
70c4762a1bSJed Brown     n    = numCells*dim;
714f99dae5SMatthew G. Knepley     if (!feFlg) {
724f99dae5SMatthew G. Knepley       ierr = PetscFree4(options->v0, options->J, options->invJ, options->detJ);CHKERRQ(ierr);
734f99dae5SMatthew G. Knepley       options->v0 = options->J = options->invJ = options->detJ = NULL;
744f99dae5SMatthew G. Knepley     }
759bf2564aSMatt McGurn     ierr = PetscOptionsRealArray("-centroid", "Input centroid for each cell", "ex8.c", options->centroid, &n, &flg);CHKERRQ(ierr);
76*98921bdaSJacob Faibussowitsch     if (flg && n != numCells*dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of centroid %D should be %D", n, numCells*dim);
779bf2564aSMatt McGurn     if (!flg) {
789bf2564aSMatt McGurn       ierr = PetscFree(options->centroid);CHKERRQ(ierr);
799bf2564aSMatt McGurn       options->centroid = NULL;
809bf2564aSMatt McGurn     }
81c4762a1bSJed Brown     n    = numCells*dim;
824f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-normal", "Input normal for each cell", "ex8.c", options->normal, &n, &flg);CHKERRQ(ierr);
83*98921bdaSJacob Faibussowitsch     if (flg && n != numCells*dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of normal %D should be %D", n, numCells*dim);
849bf2564aSMatt McGurn     if (!flg) {
859bf2564aSMatt McGurn       ierr = PetscFree(options->normal);CHKERRQ(ierr);
869bf2564aSMatt McGurn       options->normal = NULL;
879bf2564aSMatt McGurn     }
88c4762a1bSJed Brown     n    = numCells;
894f99dae5SMatthew G. Knepley     ierr = PetscOptionsRealArray("-vol", "Input volume for each cell", "ex8.c", options->vol, &n, &flg);CHKERRQ(ierr);
90*98921bdaSJacob Faibussowitsch     if (flg && n != numCells) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of vol %D should be %D", n, numCells);
919bf2564aSMatt McGurn     if (!flg) {
929bf2564aSMatt McGurn       ierr = PetscFree(options->vol);CHKERRQ(ierr);
939bf2564aSMatt McGurn       options->vol = NULL;
944f99dae5SMatthew G. Knepley     }
95c4762a1bSJed Brown   } else if (options->runType == RUN_DISPLAY) {
964f99dae5SMatthew G. Knepley     ierr = ReadMesh(PETSC_COMM_WORLD, options, &options->dm);CHKERRQ(ierr);
97c4762a1bSJed Brown   }
981e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   if (options->transform) {ierr = PetscPrintf(comm, "Using random transforms\n");CHKERRQ(ierr);}
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown static PetscErrorCode ChangeCoordinates(DM dm, PetscInt spaceDim, PetscScalar vertexCoords[])
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   PetscSection   coordSection;
107c4762a1bSJed Brown   Vec            coordinates;
108c4762a1bSJed Brown   PetscScalar   *coords;
109c4762a1bSJed Brown   PetscInt       vStart, vEnd, v, d, coordSize;
110c4762a1bSJed Brown   PetscErrorCode ierr;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   PetscFunctionBegin;
113c4762a1bSJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = PetscSectionSetChart(coordSection, vStart, vEnd);CHKERRQ(ierr);
118c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
119c4762a1bSJed Brown     ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr);
121c4762a1bSJed Brown   }
122c4762a1bSJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
129c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
130c4762a1bSJed Brown     PetscInt off;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
133c4762a1bSJed Brown     for (d = 0; d < spaceDim; ++d) {
134c4762a1bSJed Brown       coords[off+d] = vertexCoords[(v-vStart)*spaceDim+d];
135c4762a1bSJed Brown     }
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = DMSetCoordinateDim(dm, spaceDim);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
142c4762a1bSJed Brown   PetscFunctionReturn(0);
143c4762a1bSJed Brown }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown #define RelativeError(a,b) PetscAbs(a-b)/(1.0+PetscMax(PetscAbs(a),PetscAbs(b)))
146c4762a1bSJed Brown 
147c4762a1bSJed Brown static PetscErrorCode CheckFEMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx)
148c4762a1bSJed Brown {
149c4762a1bSJed Brown   PetscReal      v0[3], J[9], invJ[9], detJ;
150c4762a1bSJed Brown   PetscInt       d, i, j;
151c4762a1bSJed Brown   PetscErrorCode ierr;
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   PetscFunctionBegin;
154c4762a1bSJed Brown   ierr = DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
155c4762a1bSJed Brown   for (d = 0; d < spaceDim; ++d) {
156c4762a1bSJed Brown     if (v0[d] != v0Ex[d]) {
157c4762a1bSJed Brown       switch (spaceDim) {
158*98921bdaSJacob Faibussowitsch       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]);
159*98921bdaSJacob Faibussowitsch       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]);
160*98921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid space dimension %D", spaceDim);
161c4762a1bSJed Brown       }
162c4762a1bSJed Brown     }
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown   for (i = 0; i < spaceDim; ++i) {
165c4762a1bSJed Brown     for (j = 0; j < spaceDim; ++j) {
166*98921bdaSJacob Faibussowitsch       if (RelativeError(J[i*spaceDim+j],JEx[i*spaceDim+j])    > 10*PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid J[%D,%D]: %g != %g", i, j, (double)J[i*spaceDim+j], (double)JEx[i*spaceDim+j]);
167*98921bdaSJacob Faibussowitsch       if (RelativeError(invJ[i*spaceDim+j],invJEx[i*spaceDim+j]) > 10*PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid invJ[%D,%D]: %g != %g", i, j, (double)invJ[i*spaceDim+j], (double)invJEx[i*spaceDim+j]);
168c4762a1bSJed Brown     }
169c4762a1bSJed Brown   }
170*98921bdaSJacob Faibussowitsch   if (RelativeError(detJ,detJEx) > 10*PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid |J| = %g != %g diff %g", (double)detJ, (double)detJEx,(double)(detJ - detJEx));
171c4762a1bSJed Brown   PetscFunctionReturn(0);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown static PetscErrorCode CheckFVMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx)
175c4762a1bSJed Brown {
1764f99dae5SMatthew G. Knepley   PetscReal      tol = PetscMax(10*PETSC_SMALL, 1e-10);
177c4762a1bSJed Brown   PetscReal      centroid[3], normal[3], vol;
178c4762a1bSJed Brown   PetscInt       d;
179c4762a1bSJed Brown   PetscErrorCode ierr;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBegin;
1829bf2564aSMatt McGurn   ierr = DMPlexComputeCellGeometryFVM(dm, cell, volEx? &vol : NULL, centroidEx? centroid : NULL, normalEx? normal : NULL);CHKERRQ(ierr);
183c4762a1bSJed Brown   for (d = 0; d < spaceDim; ++d) {
1849bf2564aSMatt McGurn     if (centroidEx)
185*98921bdaSJacob Faibussowitsch       if (RelativeError(centroid[d],centroidEx[d]) > tol) SETERRQ(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]));
1869bf2564aSMatt McGurn     if (normalEx)
187*98921bdaSJacob Faibussowitsch       if (RelativeError(normal[d],normalEx[d]) > tol) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D, Invalid normal[%D]: %g != %g", cell, d, (double)normal[d], (double)normalEx[d]);
188c4762a1bSJed Brown   }
1899bf2564aSMatt McGurn   if (volEx)
190*98921bdaSJacob Faibussowitsch     if (RelativeError(volEx,vol) > tol) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D, Invalid volume = %g != %g diff %g", cell, (double)vol, (double)volEx,(double)(vol - volEx));
1914f99dae5SMatthew G. Knepley   PetscFunctionReturn(0);
1924f99dae5SMatthew G. Knepley }
1934f99dae5SMatthew G. Knepley 
1944f99dae5SMatthew G. Knepley static PetscErrorCode CheckGaussLaw(DM dm, PetscInt cell)
1954f99dae5SMatthew G. Knepley {
1964f99dae5SMatthew G. Knepley   DMPolytopeType  ct;
1974f99dae5SMatthew G. Knepley   PetscReal       tol = PetscMax(10*PETSC_SMALL, 1e-10);
1984f99dae5SMatthew G. Knepley   PetscReal       normal[3], integral[3] = {0., 0., 0.}, area;
1994f99dae5SMatthew G. Knepley   const PetscInt *cone, *ornt;
2004f99dae5SMatthew G. Knepley   PetscInt        coneSize, f, dim, cdim, d;
2014f99dae5SMatthew G. Knepley   PetscErrorCode  ierr;
2024f99dae5SMatthew G. Knepley 
2034f99dae5SMatthew G. Knepley   PetscFunctionBegin;
2044f99dae5SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2054f99dae5SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
2064f99dae5SMatthew G. Knepley   if (dim != cdim) PetscFunctionReturn(0);
2074f99dae5SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cell, &ct);CHKERRQ(ierr);
2084f99dae5SMatthew G. Knepley   if (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) PetscFunctionReturn(0);
2094f99dae5SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
2104f99dae5SMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
2114f99dae5SMatthew G. Knepley   ierr = DMPlexGetConeOrientation(dm, cell, &ornt);CHKERRQ(ierr);
2124f99dae5SMatthew G. Knepley   for (f = 0; f < coneSize; ++f) {
2139bf2564aSMatt McGurn     const PetscInt sgn = dim == 1? (f == 0 ? -1 : 1) : (ornt[f] < 0 ? -1 : 1);
2144f99dae5SMatthew G. Knepley 
2154f99dae5SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], &area, NULL, normal);CHKERRQ(ierr);
2164f99dae5SMatthew G. Knepley     for (d = 0; d < cdim; ++d) integral[d] += sgn*area*normal[d];
2174f99dae5SMatthew G. Knepley   }
218*98921bdaSJacob Faibussowitsch   for (d = 0; d < cdim; ++d) if (PetscAbsReal(integral[d]) > tol) SETERRQ(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]);
219c4762a1bSJed Brown   PetscFunctionReturn(0);
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed 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[])
223c4762a1bSJed Brown {
224c4762a1bSJed Brown   const PetscInt *cone;
225c4762a1bSJed Brown   PetscInt        coneSize, c;
226c4762a1bSJed Brown   PetscInt        dim, depth, cdim;
227c4762a1bSJed Brown   PetscErrorCode  ierr;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   PetscFunctionBegin;
230c4762a1bSJed Brown   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
233c4762a1bSJed Brown   if (v0Ex) {
234c4762a1bSJed Brown     ierr = CheckFEMGeometry(dm, cell, cdim, v0Ex, JEx, invJEx, detJEx);CHKERRQ(ierr);
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown   if (dim == depth && centroidEx) {
237c4762a1bSJed Brown     ierr = CheckFVMGeometry(dm, cell, cdim, centroidEx, normalEx, volEx);CHKERRQ(ierr);
2384f99dae5SMatthew G. Knepley     ierr = CheckGaussLaw(dm, cell);CHKERRQ(ierr);
239c4762a1bSJed Brown     if (faceCentroidEx) {
240c4762a1bSJed Brown       ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
241c4762a1bSJed Brown       ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
242c4762a1bSJed Brown       for (c = 0; c < coneSize; ++c) {
243c4762a1bSJed Brown         ierr = CheckFVMGeometry(dm, cone[c], dim, &faceCentroidEx[c*dim], &faceNormalEx[c*dim], faceVolEx[c]);CHKERRQ(ierr);
244c4762a1bSJed Brown       }
245c4762a1bSJed Brown     }
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown   if (transform) {
248c4762a1bSJed Brown     Vec          coordinates;
249c4762a1bSJed Brown     PetscSection coordSection;
250c4762a1bSJed Brown     PetscScalar *coords = NULL, *origCoords, *newCoords;
251c4762a1bSJed Brown     PetscReal   *v0ExT, *JExT, *invJExT, detJExT=0, *centroidExT, *normalExT, volExT=0;
252c4762a1bSJed Brown     PetscReal   *faceCentroidExT, *faceNormalExT, faceVolExT;
253c4762a1bSJed Brown     PetscRandom  r, ang, ang2;
254c4762a1bSJed Brown     PetscInt     coordSize, numCorners, t;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
257c4762a1bSJed Brown     ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
258c4762a1bSJed Brown     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
259c4762a1bSJed Brown     ierr = PetscMalloc2(coordSize, &origCoords, coordSize, &newCoords);CHKERRQ(ierr);
260c4762a1bSJed Brown     ierr = PetscMalloc5(cdim, &v0ExT, cdim*cdim, &JExT, cdim*cdim, &invJExT, cdim, &centroidExT, cdim, &normalExT);CHKERRQ(ierr);
261c4762a1bSJed Brown     ierr = PetscMalloc2(cdim, &faceCentroidExT, cdim, &faceNormalExT);CHKERRQ(ierr);
262c4762a1bSJed Brown     for (c = 0; c < coordSize; ++c) origCoords[c] = coords[c];
263c4762a1bSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr);
264c4762a1bSJed Brown     numCorners = coordSize/cdim;
265c4762a1bSJed Brown 
266c4762a1bSJed Brown     ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
267c4762a1bSJed Brown     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
268c4762a1bSJed Brown     ierr = PetscRandomSetInterval(r, 0.0, 10.0);CHKERRQ(ierr);
269c4762a1bSJed Brown     ierr = PetscRandomCreate(PETSC_COMM_SELF, &ang);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = PetscRandomSetFromOptions(ang);CHKERRQ(ierr);
271c4762a1bSJed Brown     ierr = PetscRandomSetInterval(ang, 0.0, 2*PETSC_PI);CHKERRQ(ierr);
272c4762a1bSJed Brown     ierr = PetscRandomCreate(PETSC_COMM_SELF, &ang2);CHKERRQ(ierr);
273c4762a1bSJed Brown     ierr = PetscRandomSetFromOptions(ang2);CHKERRQ(ierr);
274c4762a1bSJed Brown     ierr = PetscRandomSetInterval(ang2, 0.0, PETSC_PI);CHKERRQ(ierr);
275c4762a1bSJed Brown     for (t = 0; t < 1; ++t) {
276c4762a1bSJed Brown       PetscScalar trans[3];
277c4762a1bSJed Brown       PetscReal   R[9], rot[3], rotM[9];
278c4762a1bSJed Brown       PetscReal   scale, phi, theta, psi = 0.0, norm;
279c4762a1bSJed Brown       PetscInt    d, e, f, p;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown       for (c = 0; c < coordSize; ++c) newCoords[c] = origCoords[c];
282c4762a1bSJed Brown       ierr = PetscRandomGetValueReal(r, &scale);CHKERRQ(ierr);
283c4762a1bSJed Brown       ierr = PetscRandomGetValueReal(ang, &phi);CHKERRQ(ierr);
284c4762a1bSJed Brown       ierr = PetscRandomGetValueReal(ang2, &theta);CHKERRQ(ierr);
285c4762a1bSJed Brown       for (d = 0; d < cdim; ++d) {
286c4762a1bSJed Brown         ierr = PetscRandomGetValue(r, &trans[d]);CHKERRQ(ierr);
287c4762a1bSJed Brown       }
288c4762a1bSJed Brown       switch (cdim) {
289c4762a1bSJed Brown       case 2:
290c4762a1bSJed Brown         R[0] = PetscCosReal(phi); R[1] = -PetscSinReal(phi);
291c4762a1bSJed Brown         R[2] = PetscSinReal(phi); R[3] =  PetscCosReal(phi);
292c4762a1bSJed Brown         break;
293c4762a1bSJed Brown       case 3:
294c4762a1bSJed Brown       {
295c4762a1bSJed Brown         const PetscReal ct = PetscCosReal(theta), st = PetscSinReal(theta);
296c4762a1bSJed Brown         const PetscReal cp = PetscCosReal(phi),   sp = PetscSinReal(phi);
297c4762a1bSJed Brown         const PetscReal cs = PetscCosReal(psi),   ss = PetscSinReal(psi);
298c4762a1bSJed Brown         R[0] = ct*cs; R[1] = sp*st*cs - cp*ss;    R[2] = sp*ss    + cp*st*cs;
299c4762a1bSJed Brown         R[3] = ct*ss; R[4] = cp*cs    + sp*st*ss; R[5] = cp*st*ss - sp*cs;
300c4762a1bSJed Brown         R[6] = -st;   R[7] = sp*ct;               R[8] = cp*ct;
301c4762a1bSJed Brown         break;
302c4762a1bSJed Brown       }
303*98921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid coordinate dimension %D", cdim);
304c4762a1bSJed Brown       }
305c4762a1bSJed Brown       if (v0Ex) {
306c4762a1bSJed Brown         detJExT = detJEx;
307c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
308c4762a1bSJed Brown           v0ExT[d] = v0Ex[d];
309c4762a1bSJed Brown           for (e = 0; e < cdim; ++e) {
310c4762a1bSJed Brown             JExT[d*cdim+e]    = JEx[d*cdim+e];
311c4762a1bSJed Brown             invJExT[d*cdim+e] = invJEx[d*cdim+e];
312c4762a1bSJed Brown           }
313c4762a1bSJed Brown         }
314c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
315c4762a1bSJed Brown           v0ExT[d] *= scale;
316c4762a1bSJed Brown           v0ExT[d] += PetscRealPart(trans[d]);
317c4762a1bSJed Brown           /* Only scale dimensions in the manifold */
318c4762a1bSJed Brown           for (e = 0; e < dim; ++e) {
319c4762a1bSJed Brown             JExT[d*cdim+e]    *= scale;
320c4762a1bSJed Brown             invJExT[d*cdim+e] /= scale;
321c4762a1bSJed Brown           }
322c4762a1bSJed Brown           if (d < dim) detJExT *= scale;
323c4762a1bSJed Brown         }
324c4762a1bSJed Brown         /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */
325c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
326c4762a1bSJed Brown           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
327c4762a1bSJed Brown             rot[d] += R[d*cdim+e] * v0ExT[e];
328c4762a1bSJed Brown           }
329c4762a1bSJed Brown         }
330c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) v0ExT[d] = rot[d];
331c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
332c4762a1bSJed Brown           for (e = 0; e < cdim; ++e) {
333c4762a1bSJed Brown             for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) {
334c4762a1bSJed Brown               rotM[d*cdim+e] += R[d*cdim+f] * JExT[f*cdim+e];
335c4762a1bSJed Brown             }
336c4762a1bSJed Brown           }
337c4762a1bSJed Brown         }
338c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
339c4762a1bSJed Brown           for (e = 0; e < cdim; ++e) {
340c4762a1bSJed Brown             JExT[d*cdim+e] = rotM[d*cdim+e];
341c4762a1bSJed Brown           }
342c4762a1bSJed Brown         }
343c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
344c4762a1bSJed Brown           for (e = 0; e < cdim; ++e) {
345c4762a1bSJed Brown             for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) {
346c4762a1bSJed Brown               rotM[d*cdim+e] += invJExT[d*cdim+f] * R[e*cdim+f];
347c4762a1bSJed Brown             }
348c4762a1bSJed Brown           }
349c4762a1bSJed Brown         }
350c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
351c4762a1bSJed Brown           for (e = 0; e < cdim; ++e) {
352c4762a1bSJed Brown             invJExT[d*cdim+e] = rotM[d*cdim+e];
353c4762a1bSJed Brown           }
354c4762a1bSJed Brown         }
355c4762a1bSJed Brown       }
356c4762a1bSJed Brown       if (centroidEx) {
357c4762a1bSJed Brown         volExT = volEx;
358c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
359c4762a1bSJed Brown           centroidExT[d]  = centroidEx[d];
360c4762a1bSJed Brown           normalExT[d]    = normalEx[d];
361c4762a1bSJed Brown         }
362c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
363c4762a1bSJed Brown           centroidExT[d] *= scale;
364c4762a1bSJed Brown           centroidExT[d] += PetscRealPart(trans[d]);
365c4762a1bSJed Brown           normalExT[d]   /= scale;
366c4762a1bSJed Brown           /* Only scale dimensions in the manifold */
367c4762a1bSJed Brown           if (d < dim) volExT  *= scale;
368c4762a1bSJed Brown         }
369c4762a1bSJed Brown         /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */
370c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
371c4762a1bSJed Brown           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
372c4762a1bSJed Brown             rot[d] += R[d*cdim+e] * centroidExT[e];
373c4762a1bSJed Brown           }
374c4762a1bSJed Brown         }
375c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) centroidExT[d] = rot[d];
376c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
377c4762a1bSJed Brown           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
378c4762a1bSJed Brown             rot[d] += R[d*cdim+e] * normalExT[e];
379c4762a1bSJed Brown           }
380c4762a1bSJed Brown         }
381c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) normalExT[d] = rot[d];
382c4762a1bSJed Brown         for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(normalExT[d]);
383c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
384c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) normalExT[d] /= norm;
385c4762a1bSJed Brown       }
386c4762a1bSJed Brown       for (d = 0; d < cdim; ++d) {
387c4762a1bSJed Brown         for (p = 0; p < numCorners; ++p) {
388c4762a1bSJed Brown           newCoords[p*cdim+d] *= scale;
389c4762a1bSJed Brown           newCoords[p*cdim+d] += trans[d];
390c4762a1bSJed Brown         }
391c4762a1bSJed Brown       }
392c4762a1bSJed Brown       for (p = 0; p < numCorners; ++p) {
393c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) {
394c4762a1bSJed Brown           for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
395c4762a1bSJed Brown             rot[d] += R[d*cdim+e] * PetscRealPart(newCoords[p*cdim+e]);
396c4762a1bSJed Brown           }
397c4762a1bSJed Brown         }
398c4762a1bSJed Brown         for (d = 0; d < cdim; ++d) newCoords[p*cdim+d] = rot[d];
399c4762a1bSJed Brown       }
400c4762a1bSJed Brown 
401c4762a1bSJed Brown       ierr = ChangeCoordinates(dm, cdim, newCoords);CHKERRQ(ierr);
402c4762a1bSJed Brown       if (v0Ex) {
403c4762a1bSJed Brown         ierr = CheckFEMGeometry(dm, 0, cdim, v0ExT, JExT, invJExT, detJExT);CHKERRQ(ierr);
404c4762a1bSJed Brown       }
405c4762a1bSJed Brown       if (dim == depth && centroidEx) {
406c4762a1bSJed Brown         ierr = CheckFVMGeometry(dm, cell, cdim, centroidExT, normalExT, volExT);CHKERRQ(ierr);
4074f99dae5SMatthew G. Knepley         ierr = CheckGaussLaw(dm, cell);CHKERRQ(ierr);
408c4762a1bSJed Brown         if (faceCentroidEx) {
409c4762a1bSJed Brown           ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
410c4762a1bSJed Brown           ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
411c4762a1bSJed Brown           for (c = 0; c < coneSize; ++c) {
412c4762a1bSJed Brown             PetscInt off = c*cdim;
413c4762a1bSJed Brown 
414c4762a1bSJed Brown             faceVolExT = faceVolEx[c];
415c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) {
416c4762a1bSJed Brown               faceCentroidExT[d]  = faceCentroidEx[off+d];
417c4762a1bSJed Brown               faceNormalExT[d]    = faceNormalEx[off+d];
418c4762a1bSJed Brown             }
419c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) {
420c4762a1bSJed Brown               faceCentroidExT[d] *= scale;
421c4762a1bSJed Brown               faceCentroidExT[d] += PetscRealPart(trans[d]);
422c4762a1bSJed Brown               faceNormalExT[d]   /= scale;
423c4762a1bSJed Brown               /* Only scale dimensions in the manifold */
424c4762a1bSJed Brown               if (d < dim-1) {
425c4762a1bSJed Brown                 faceVolExT *= scale;
426c4762a1bSJed Brown               }
427c4762a1bSJed Brown             }
428c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) {
429c4762a1bSJed Brown               for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
430c4762a1bSJed Brown                 rot[d] += R[d*cdim+e] * faceCentroidExT[e];
431c4762a1bSJed Brown               }
432c4762a1bSJed Brown             }
433c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) faceCentroidExT[d] = rot[d];
434c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) {
435c4762a1bSJed Brown               for (e = 0, rot[d] = 0.0; e < cdim; ++e) {
436c4762a1bSJed Brown                 rot[d] += R[d*cdim+e] * faceNormalExT[e];
437c4762a1bSJed Brown               }
438c4762a1bSJed Brown             }
439c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) faceNormalExT[d] = rot[d];
440c4762a1bSJed Brown             for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(faceNormalExT[d]);
441c4762a1bSJed Brown             norm = PetscSqrtReal(norm);
442c4762a1bSJed Brown             for (d = 0; d < cdim; ++d) faceNormalExT[d] /= norm;
443c4762a1bSJed Brown 
444c4762a1bSJed Brown             ierr = CheckFVMGeometry(dm, cone[c], cdim, faceCentroidExT, faceNormalExT, faceVolExT);CHKERRQ(ierr);
445c4762a1bSJed Brown           }
446c4762a1bSJed Brown         }
447c4762a1bSJed Brown       }
448c4762a1bSJed Brown     }
449c4762a1bSJed Brown     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
450c4762a1bSJed Brown     ierr = PetscRandomDestroy(&ang);CHKERRQ(ierr);
451c4762a1bSJed Brown     ierr = PetscRandomDestroy(&ang2);CHKERRQ(ierr);
452c4762a1bSJed Brown     ierr = PetscFree2(origCoords, newCoords);CHKERRQ(ierr);
453c4762a1bSJed Brown     ierr = PetscFree5(v0ExT, JExT, invJExT, centroidExT, normalExT);CHKERRQ(ierr);
454c4762a1bSJed Brown     ierr = PetscFree2(faceCentroidExT, faceNormalExT);CHKERRQ(ierr);
455c4762a1bSJed Brown   }
456c4762a1bSJed Brown   PetscFunctionReturn(0);
457c4762a1bSJed Brown }
458c4762a1bSJed Brown 
4594f99dae5SMatthew G. Knepley static PetscErrorCode TestTriangle(MPI_Comm comm, PetscBool transform)
460c4762a1bSJed Brown {
461c4762a1bSJed Brown   DM             dm;
462c4762a1bSJed Brown   PetscErrorCode ierr;
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   PetscFunctionBegin;
4654f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRIANGLE, &dm);CHKERRQ(ierr);
466c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
467c4762a1bSJed Brown   /* Check reference geometry: determinant is scaled by reference volume (2.0) */
468c4762a1bSJed Brown   {
469c4762a1bSJed Brown     PetscReal v0Ex[2]       = {-1.0, -1.0};
470c4762a1bSJed Brown     PetscReal JEx[4]        = {1.0, 0.0, 0.0, 1.0};
471c4762a1bSJed Brown     PetscReal invJEx[4]     = {1.0, 0.0, 0.0, 1.0};
472c4762a1bSJed Brown     PetscReal detJEx        = 1.0;
473c4762a1bSJed Brown     PetscReal centroidEx[2] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.)};
474c4762a1bSJed Brown     PetscReal normalEx[2]   = {0.0, 0.0};
475c4762a1bSJed Brown     PetscReal volEx         = 2.0;
476c4762a1bSJed Brown 
477c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
478c4762a1bSJed Brown   }
479c4762a1bSJed Brown   /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (2.0) */
480c4762a1bSJed Brown   {
481c4762a1bSJed Brown     PetscScalar vertexCoords[9] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, -1.0, 1.0, 0.0};
482c4762a1bSJed Brown     PetscReal   v0Ex[3]         = {-1.0, -1.0, 0.0};
483c4762a1bSJed Brown     PetscReal   JEx[9]          = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
484c4762a1bSJed Brown     PetscReal   invJEx[9]       = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
485c4762a1bSJed Brown     PetscReal   detJEx          = 1.0;
486c4762a1bSJed Brown     PetscReal   centroidEx[3]   = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0};
487c4762a1bSJed Brown     PetscReal   normalEx[3]     = {0.0, 0.0, 1.0};
488c4762a1bSJed Brown     PetscReal   volEx           = 2.0;
489c4762a1bSJed Brown 
4904f99dae5SMatthew G. Knepley     ierr = ChangeCoordinates(dm, 3, vertexCoords);CHKERRQ(ierr);
491c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
492c4762a1bSJed Brown   }
493c4762a1bSJed Brown   /* Cleanup */
494c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
495c4762a1bSJed Brown   PetscFunctionReturn(0);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
4984f99dae5SMatthew G. Knepley static PetscErrorCode TestQuadrilateral(MPI_Comm comm, PetscBool transform)
499c4762a1bSJed Brown {
500c4762a1bSJed Brown   DM             dm;
501c4762a1bSJed Brown   PetscErrorCode ierr;
502c4762a1bSJed Brown 
503c4762a1bSJed Brown   PetscFunctionBegin;
5044f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_QUADRILATERAL, &dm);CHKERRQ(ierr);
505c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
506c4762a1bSJed Brown   /* Check reference geometry: determinant is scaled by reference volume (2.0) */
507c4762a1bSJed Brown   {
508c4762a1bSJed Brown     PetscReal v0Ex[2]       = {-1.0, -1.0};
509c4762a1bSJed Brown     PetscReal JEx[4]        = {1.0, 0.0, 0.0, 1.0};
510c4762a1bSJed Brown     PetscReal invJEx[4]     = {1.0, 0.0, 0.0, 1.0};
511c4762a1bSJed Brown     PetscReal detJEx        = 1.0;
512c4762a1bSJed Brown     PetscReal centroidEx[2] = {0.0, 0.0};
513c4762a1bSJed Brown     PetscReal normalEx[2]   = {0.0, 0.0};
514c4762a1bSJed Brown     PetscReal volEx         = 4.0;
515c4762a1bSJed Brown 
516c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
517c4762a1bSJed Brown   }
518c4762a1bSJed Brown   /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (4.0) */
519c4762a1bSJed Brown   {
520c4762a1bSJed 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};
521c4762a1bSJed Brown     PetscReal   v0Ex[3]          = {-1.0, -1.0, 0.0};
522c4762a1bSJed Brown     PetscReal   JEx[9]           = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
523c4762a1bSJed Brown     PetscReal   invJEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
524c4762a1bSJed Brown     PetscReal   detJEx           = 1.0;
525c4762a1bSJed Brown     PetscReal   centroidEx[3]    = {0.0, 0.0, 0.0};
526c4762a1bSJed Brown     PetscReal   normalEx[3]      = {0.0, 0.0, 1.0};
527c4762a1bSJed Brown     PetscReal   volEx            = 4.0;
528c4762a1bSJed Brown 
5294f99dae5SMatthew G. Knepley     ierr = ChangeCoordinates(dm, 3, vertexCoords);CHKERRQ(ierr);
530c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
531c4762a1bSJed Brown   }
532c4762a1bSJed Brown   /* Cleanup */
533c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
534c4762a1bSJed Brown   PetscFunctionReturn(0);
535c4762a1bSJed Brown }
536c4762a1bSJed Brown 
5374f99dae5SMatthew G. Knepley static PetscErrorCode TestTetrahedron(MPI_Comm comm, PetscBool transform)
538c4762a1bSJed Brown {
539c4762a1bSJed Brown   DM             dm;
540c4762a1bSJed Brown   PetscErrorCode ierr;
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   PetscFunctionBegin;
5434f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TETRAHEDRON, &dm);CHKERRQ(ierr);
544c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
545c4762a1bSJed Brown   /* Check reference geometry: determinant is scaled by reference volume (4/3) */
546c4762a1bSJed Brown   {
547c4762a1bSJed Brown     PetscReal v0Ex[3]       = {-1.0, -1.0, -1.0};
548c4762a1bSJed Brown     PetscReal JEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
549c4762a1bSJed Brown     PetscReal invJEx[9]     = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
550c4762a1bSJed Brown     PetscReal detJEx        = 1.0;
551c4762a1bSJed Brown     PetscReal centroidEx[3] = {-0.5, -0.5, -0.5};
552c4762a1bSJed Brown     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
553c4762a1bSJed Brown     PetscReal volEx         = (PetscReal)4.0/(PetscReal)3.0;
554c4762a1bSJed Brown 
555c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
556c4762a1bSJed Brown   }
557c4762a1bSJed Brown   /* Cleanup */
558c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
559c4762a1bSJed Brown   PetscFunctionReturn(0);
560c4762a1bSJed Brown }
561c4762a1bSJed Brown 
5624f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedron(MPI_Comm comm, PetscBool transform)
563c4762a1bSJed Brown {
564c4762a1bSJed Brown   DM             dm;
565c4762a1bSJed Brown   PetscErrorCode ierr;
566c4762a1bSJed Brown 
567c4762a1bSJed Brown   PetscFunctionBegin;
5684f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm);CHKERRQ(ierr);
569c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
570c4762a1bSJed Brown   /* Check reference geometry: determinant is scaled by reference volume 8.0 */
571c4762a1bSJed Brown   {
572c4762a1bSJed Brown     PetscReal v0Ex[3]       = {-1.0, -1.0, -1.0};
573c4762a1bSJed Brown     PetscReal JEx[9]        = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
574c4762a1bSJed Brown     PetscReal invJEx[9]     = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
575c4762a1bSJed Brown     PetscReal detJEx        = 1.0;
576c4762a1bSJed Brown     PetscReal centroidEx[3] = {0.0, 0.0, 0.0};
577c4762a1bSJed Brown     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
578c4762a1bSJed Brown     PetscReal volEx         = 8.0;
579c4762a1bSJed Brown 
580c4762a1bSJed Brown     ierr = CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
581c4762a1bSJed Brown   }
582c4762a1bSJed Brown   /* Cleanup */
583c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
584c4762a1bSJed Brown   PetscFunctionReturn(0);
585c4762a1bSJed Brown }
586c4762a1bSJed Brown 
5874f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedronCurved(MPI_Comm comm)
588c4762a1bSJed Brown {
589c4762a1bSJed Brown   DM             dm;
5904f99dae5SMatthew 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,
5914f99dae5SMatthew 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};
592c4762a1bSJed Brown   PetscErrorCode ierr;
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   PetscFunctionBegin;
5954f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm);CHKERRQ(ierr);
5964f99dae5SMatthew G. Knepley   ierr = ChangeCoordinates(dm, 3, coords);CHKERRQ(ierr);
597c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
5984f99dae5SMatthew G. Knepley   {
5994f99dae5SMatthew G. Knepley     PetscReal centroidEx[3] = {0.0, 0.0, 0.016803278688524603};
6004f99dae5SMatthew G. Knepley     PetscReal normalEx[3]   = {0.0, 0.0, 0.0};
6014f99dae5SMatthew G. Knepley     PetscReal volEx         = 8.1333333333333346;
6024f99dae5SMatthew G. Knepley 
6034f99dae5SMatthew G. Knepley     ierr = CheckCell(dm, 0, PETSC_FALSE, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, NULL, NULL, NULL);CHKERRQ(ierr);
604c4762a1bSJed Brown   }
6054f99dae5SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
6064f99dae5SMatthew G. Knepley   PetscFunctionReturn(0);
6074f99dae5SMatthew G. Knepley }
6084f99dae5SMatthew G. Knepley 
6094f99dae5SMatthew G. Knepley /* This wedge is a tensor product cell, rather than a normal wedge */
6104f99dae5SMatthew G. Knepley static PetscErrorCode TestWedge(MPI_Comm comm, PetscBool transform)
6114f99dae5SMatthew G. Knepley {
6124f99dae5SMatthew G. Knepley   DM             dm;
6134f99dae5SMatthew G. Knepley   PetscErrorCode ierr;
6144f99dae5SMatthew G. Knepley 
6154f99dae5SMatthew G. Knepley   PetscFunctionBegin;
6164f99dae5SMatthew G. Knepley   ierr = DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRI_PRISM_TENSOR, &dm);CHKERRQ(ierr);
6174f99dae5SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
618c4762a1bSJed Brown   /* Check reference geometry: determinant is scaled by reference volume 4.0 */
619c4762a1bSJed Brown   {
620c4762a1bSJed Brown #if 0
621c4762a1bSJed Brown     /* FEM geometry is not functional for wedges */
622c4762a1bSJed Brown     PetscReal v0Ex[3]   = {-1.0, -1.0, -1.0};
623c4762a1bSJed Brown     PetscReal JEx[9]    = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
624c4762a1bSJed Brown     PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0};
625c4762a1bSJed Brown     PetscReal detJEx    = 1.0;
626c4762a1bSJed Brown #endif
627c4762a1bSJed Brown 
6284f99dae5SMatthew G. Knepley     {
629c4762a1bSJed Brown       PetscReal       centroidEx[3]      = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0};
630c4762a1bSJed Brown       PetscReal       normalEx[3]        = {0.0, 0.0, 0.0};
631c4762a1bSJed Brown       PetscReal       volEx              = 4.0;
632c4762a1bSJed Brown       PetscReal       faceVolEx[5]       = {2.0, 2.0, 4.0, PETSC_SQRT2*4.0, 4.0};
633c4762a1bSJed 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};
634c4762a1bSJed Brown       PetscReal       faceCentroidEx[15] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), -1.0,
635c4762a1bSJed Brown                                             -((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.),  1.0,
636c4762a1bSJed Brown                                             0.0, -1.0, 0.0,  0.0, 0.0, 0.0,  -1.0, 0.0, 0.0};
637c4762a1bSJed Brown 
638c4762a1bSJed Brown       ierr = CheckCell(dm, 0, transform, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, faceCentroidEx, faceNormalEx, faceVolEx);CHKERRQ(ierr);
639c4762a1bSJed Brown     }
640c4762a1bSJed Brown   }
641c4762a1bSJed Brown   /* Cleanup */
642c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
643c4762a1bSJed Brown   PetscFunctionReturn(0);
644c4762a1bSJed Brown }
645c4762a1bSJed Brown 
646c4762a1bSJed Brown int main(int argc, char **argv)
647c4762a1bSJed Brown {
648c4762a1bSJed Brown   AppCtx         user;
649c4762a1bSJed Brown   PetscErrorCode ierr;
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
652c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
653c4762a1bSJed Brown   if (user.runType == RUN_REFERENCE) {
6544f99dae5SMatthew G. Knepley     ierr = TestTriangle(PETSC_COMM_SELF, user.transform);CHKERRQ(ierr);
6554f99dae5SMatthew G. Knepley     ierr = TestQuadrilateral(PETSC_COMM_SELF, user.transform);CHKERRQ(ierr);
6564f99dae5SMatthew G. Knepley     ierr = TestTetrahedron(PETSC_COMM_SELF, user.transform);CHKERRQ(ierr);
6574f99dae5SMatthew G. Knepley     ierr = TestHexahedron(PETSC_COMM_SELF, user.transform);CHKERRQ(ierr);
6584f99dae5SMatthew G. Knepley     ierr = TestWedge(PETSC_COMM_SELF, user.transform);CHKERRQ(ierr);
6594f99dae5SMatthew G. Knepley   } else if (user.runType == RUN_HEX_CURVED) {
6604f99dae5SMatthew G. Knepley     ierr = TestHexahedronCurved(PETSC_COMM_SELF);CHKERRQ(ierr);
661c4762a1bSJed Brown   } else if (user.runType == RUN_FILE) {
662c4762a1bSJed Brown     PetscInt dim, cStart, cEnd, c;
663c4762a1bSJed Brown 
664c4762a1bSJed Brown     ierr = DMGetDimension(user.dm, &dim);CHKERRQ(ierr);
665c4762a1bSJed Brown     ierr = DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
666c4762a1bSJed Brown     for (c = 0; c < cEnd-cStart; ++c) {
6674f99dae5SMatthew G. Knepley       PetscReal *v0       = user.v0       ? &user.v0[c*dim] : NULL;
6684f99dae5SMatthew G. Knepley       PetscReal *J        = user.J        ? &user.J[c*dim*dim] : NULL;
6694f99dae5SMatthew G. Knepley       PetscReal *invJ     = user.invJ     ? &user.invJ[c*dim*dim] : NULL;
6704f99dae5SMatthew G. Knepley       PetscReal  detJ     = user.detJ     ?  user.detJ[c] : 0.0;
6714f99dae5SMatthew G. Knepley       PetscReal *centroid = user.centroid ? &user.centroid[c*dim] : NULL;
6724f99dae5SMatthew G. Knepley       PetscReal *normal   = user.normal   ? &user.normal[c*dim] : NULL;
6734f99dae5SMatthew G. Knepley       PetscReal  vol      = user.vol      ?  user.vol[c] : 0.0;
6744f99dae5SMatthew G. Knepley 
6754f99dae5SMatthew G. Knepley       ierr = CheckCell(user.dm, c+cStart, PETSC_FALSE, v0, J, invJ, detJ, centroid, normal, vol, NULL, NULL, NULL);CHKERRQ(ierr);
676c4762a1bSJed Brown     }
6774f99dae5SMatthew G. Knepley     ierr = PetscFree4(user.v0,user.J,user.invJ,user.detJ);CHKERRQ(ierr);
6789bf2564aSMatt McGurn     ierr = PetscFree(user.centroid);CHKERRQ(ierr);
6799bf2564aSMatt McGurn     ierr = PetscFree(user.normal);CHKERRQ(ierr);
6809bf2564aSMatt McGurn     ierr = PetscFree(user.vol);CHKERRQ(ierr);
681c4762a1bSJed Brown     ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
682c4762a1bSJed Brown   } else if (user.runType == RUN_DISPLAY) {
683c4762a1bSJed Brown     DM                 gdm, dmCell;
684c4762a1bSJed Brown     Vec                cellgeom, facegeom;
685c4762a1bSJed Brown     const PetscScalar *cgeom;
686c4762a1bSJed Brown     PetscInt           dim, d, cStart, cEnd, cEndInterior, c;
687c4762a1bSJed Brown 
688c4762a1bSJed Brown     ierr = DMGetCoordinateDim(user.dm, &dim);CHKERRQ(ierr);
689c4762a1bSJed Brown     ierr = DMPlexConstructGhostCells(user.dm, NULL, NULL, &gdm);CHKERRQ(ierr);
690c4762a1bSJed Brown     if (gdm) {
691c4762a1bSJed Brown       ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
692c4762a1bSJed Brown       user.dm = gdm;
693c4762a1bSJed Brown     }
694c4762a1bSJed Brown     ierr = DMPlexComputeGeometryFVM(user.dm, &cellgeom, &facegeom);CHKERRQ(ierr);
695c4762a1bSJed Brown     ierr = DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
696c4762a1bSJed Brown     ierr = DMPlexGetGhostCellStratum(user.dm, &cEndInterior, NULL);CHKERRQ(ierr);
697c4762a1bSJed Brown     if (cEndInterior >= 0) cEnd = cEndInterior;
698c4762a1bSJed Brown     ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr);
699c4762a1bSJed Brown     ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr);
700c4762a1bSJed Brown     for (c = 0; c < cEnd-cStart; ++c) {
701c4762a1bSJed Brown       PetscFVCellGeom *cg;
702c4762a1bSJed Brown 
703c4762a1bSJed Brown       ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
704c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %4D: Centroid (", c);CHKERRQ(ierr);
705c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
706c4762a1bSJed Brown         if (d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);}
707c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "%12.2g", cg->centroid[d]);CHKERRQ(ierr);
708c4762a1bSJed Brown       }
709c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, ") Vol %12.2g\n", cg->volume);CHKERRQ(ierr);
710c4762a1bSJed Brown     }
711c4762a1bSJed Brown     ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr);
712c4762a1bSJed Brown     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
713c4762a1bSJed Brown     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
714c4762a1bSJed Brown     ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
715c4762a1bSJed Brown   }
716c4762a1bSJed Brown   ierr = PetscFinalize();
717c4762a1bSJed Brown   return ierr;
718c4762a1bSJed Brown }
719c4762a1bSJed Brown 
720c4762a1bSJed Brown /*TEST
721c4762a1bSJed Brown 
722c4762a1bSJed Brown   test:
7234f99dae5SMatthew G. Knepley     suffix: 1
724c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail
725c4762a1bSJed Brown   test:
726c4762a1bSJed Brown     suffix: 2
7274f99dae5SMatthew G. Knepley     args: -run_type hex_curved
728c4762a1bSJed Brown   test:
729c4762a1bSJed Brown     suffix: 3
7304f99dae5SMatthew G. Knepley     args: -transform
731c4762a1bSJed Brown   test:
732c4762a1bSJed Brown     suffix: 4
733c4762a1bSJed Brown     requires: exodusii
7344f99dae5SMatthew 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
735c4762a1bSJed Brown   test:
736c4762a1bSJed Brown     suffix: 5
7374f99dae5SMatthew 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
7389bf2564aSMatt McGurn   test:
7399bf2564aSMatt McGurn     suffix: 6
7409bf2564aSMatt McGurn     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
741c4762a1bSJed Brown TEST*/
742