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