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