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