18b438e21SMatthew G. Knepley static char help[] = "Tests for high order geometry\n\n"; 28b438e21SMatthew G. Knepley 38b438e21SMatthew G. Knepley #include <petscdmplex.h> 48b438e21SMatthew G. Knepley #include <petscds.h> 58b438e21SMatthew G. Knepley 660c1a66aSMatthew G. Knepley typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_FLARE, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform; 760c1a66aSMatthew G. Knepley const char * const TransformTypes[] = {"none", "shear", "flare", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL}; 88b438e21SMatthew G. Knepley 98b438e21SMatthew G. Knepley typedef struct { 1030602db0SMatthew G. Knepley PetscBool coordSpace; /* Flag to create coordinate space */ 118b438e21SMatthew G. Knepley Transform meshTransform; /* Transform for initial box mesh */ 128b438e21SMatthew G. Knepley PetscReal *transformDataReal; /* Parameters for mesh transform */ 138b438e21SMatthew G. Knepley PetscScalar *transformData; /* Parameters for mesh transform */ 148b438e21SMatthew G. Knepley PetscReal volume; /* Analytical volume of the mesh */ 158b438e21SMatthew G. Knepley PetscReal tol; /* Tolerance for volume check */ 168b438e21SMatthew G. Knepley } AppCtx; 178b438e21SMatthew G. Knepley 188b438e21SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 198b438e21SMatthew G. Knepley { 208b438e21SMatthew G. Knepley PetscInt n = 0, i; 218b438e21SMatthew G. Knepley PetscErrorCode ierr; 228b438e21SMatthew G. Knepley 238b438e21SMatthew G. Knepley PetscFunctionBegin; 2430602db0SMatthew G. Knepley options->coordSpace = PETSC_TRUE; 258b438e21SMatthew G. Knepley options->meshTransform = TRANSFORM_NONE; 268b438e21SMatthew G. Knepley options->transformDataReal = NULL; 278b438e21SMatthew G. Knepley options->transformData = NULL; 288b438e21SMatthew G. Knepley options->volume = -1.0; 298b438e21SMatthew G. Knepley options->tol = PETSC_SMALL; 308b438e21SMatthew G. Knepley 318b438e21SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr); 3230602db0SMatthew G. Knepley ierr = PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL);CHKERRQ(ierr); 338b438e21SMatthew G. Knepley ierr = PetscOptionsEnum("-mesh_transform", "Method to transform initial box mesh <none, shear, annulus, shell>", "ex33.c", TransformTypes, (PetscEnum) options->meshTransform, (PetscEnum *) &options->meshTransform, NULL);CHKERRQ(ierr); 348b438e21SMatthew G. Knepley switch (options->meshTransform) { 358b438e21SMatthew G. Knepley case TRANSFORM_NONE: break; 368b438e21SMatthew G. Knepley case TRANSFORM_SHEAR: 3739290ff6SMatthew G. Knepley n = 2; 388b438e21SMatthew G. Knepley ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr); 398b438e21SMatthew G. Knepley for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0; 408b438e21SMatthew G. Knepley ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr); 418b438e21SMatthew G. Knepley break; 4260c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 4360c1a66aSMatthew G. Knepley n = 4; 4460c1a66aSMatthew G. Knepley ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr); 4560c1a66aSMatthew G. Knepley for (i = 0; i < n; ++i) options->transformData[i] = 1.0; 4660c1a66aSMatthew G. Knepley options->transformData[0] = (PetscScalar) 0; 4760c1a66aSMatthew G. Knepley ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr); 4860c1a66aSMatthew G. Knepley break; 498b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 508b438e21SMatthew G. Knepley n = 2; 518b438e21SMatthew G. Knepley ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr); 528b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 538b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 548b438e21SMatthew G. Knepley ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr); 558b438e21SMatthew G. Knepley break; 568b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 578b438e21SMatthew G. Knepley n = 2; 588b438e21SMatthew G. Knepley ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr); 598b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 608b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 618b438e21SMatthew G. Knepley ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr); 628b438e21SMatthew G. Knepley break; 63*98921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform); 648b438e21SMatthew G. Knepley } 658b438e21SMatthew G. Knepley ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr); 668b438e21SMatthew G. Knepley ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr); 671e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 688b438e21SMatthew G. Knepley PetscFunctionReturn(0); 698b438e21SMatthew G. Knepley } 708b438e21SMatthew G. Knepley 718b438e21SMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 728b438e21SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 738b438e21SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 748b438e21SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 758b438e21SMatthew G. Knepley { 768b438e21SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 778b438e21SMatthew G. Knepley PetscInt c; 788b438e21SMatthew G. Knepley 798b438e21SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 808b438e21SMatthew G. Knepley } 818b438e21SMatthew G. Knepley 8260c1a66aSMatthew G. Knepley /* Flare applies the transformation, assuming we fix x_f, 8360c1a66aSMatthew G. Knepley 8460c1a66aSMatthew G. Knepley x_i = x_i * alpha_i x_f 8560c1a66aSMatthew G. Knepley */ 8660c1a66aSMatthew G. Knepley static void f0_flare(PetscInt dim, PetscInt Nf, PetscInt NfAux, 8760c1a66aSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 8860c1a66aSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 8960c1a66aSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[]) 9060c1a66aSMatthew G. Knepley { 9160c1a66aSMatthew G. Knepley const PetscInt Nc = uOff[1]-uOff[0]; 9260c1a66aSMatthew G. Knepley const PetscInt cf = (PetscInt) PetscRealPart(constants[0]); 9360c1a66aSMatthew G. Knepley PetscInt c; 9460c1a66aSMatthew G. Knepley 9560c1a66aSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 9660c1a66aSMatthew G. Knepley coords[c] = u[c] * (c == cf ? 1.0 : constants[c+1]*u[cf]); 9760c1a66aSMatthew G. Knepley } 9860c1a66aSMatthew G. Knepley } 9960c1a66aSMatthew G. Knepley 1008b438e21SMatthew G. Knepley /* 1018b438e21SMatthew G. Knepley We would like to map the unit square to a quarter of the annulus between circles of radius 1 and 2. We start by mapping the straight sections, which 1028b438e21SMatthew G. Knepley will correspond to the top and bottom of our square. So 1038b438e21SMatthew G. Knepley 1048b438e21SMatthew G. Knepley (0,0)--(1,0) ==> (1,0)--(2,0) Just a shift of (1,0) 1058b438e21SMatthew G. Knepley (0,1)--(1,1) ==> (0,1)--(0,2) Switch x and y 1068b438e21SMatthew G. Knepley 1078b438e21SMatthew G. Knepley So it looks like we want to map each layer in y to a ray, so x is the radius and y is the angle: 1088b438e21SMatthew G. Knepley 1098b438e21SMatthew G. Knepley (x, y) ==> (x+1, \pi/2 y) in (r', \theta') space 1108b438e21SMatthew G. Knepley ==> ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space 1118b438e21SMatthew G. Knepley */ 1128b438e21SMatthew G. Knepley static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1138b438e21SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1148b438e21SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1158b438e21SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 1168b438e21SMatthew G. Knepley { 1178b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1188b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1198b438e21SMatthew G. Knepley 1208b438e21SMatthew G. Knepley xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]); 1218b438e21SMatthew G. Knepley xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]); 1228b438e21SMatthew G. Knepley } 1238b438e21SMatthew G. Knepley 1248b438e21SMatthew G. Knepley /* 1258b438e21SMatthew G. Knepley We would like to map the unit cube to a hemisphere of the spherical shell between balls of radius 1 and 2. We want to map the bottom surface onto the 1268b438e21SMatthew G. Knepley lower hemisphere and the upper surface onto the top, letting z be the radius. 1278b438e21SMatthew G. Knepley 1288b438e21SMatthew G. Knepley (x, y) ==> ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x) in (r', \theta', \phi') space 1298b438e21SMatthew G. Knepley ==> ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space 1308b438e21SMatthew G. Knepley */ 1318b438e21SMatthew G. Knepley static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1328b438e21SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1338b438e21SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1348b438e21SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 1358b438e21SMatthew G. Knepley { 1368b438e21SMatthew G. Knepley const PetscReal pi4 = PETSC_PI/4.0; 1378b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1388b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1398b438e21SMatthew G. Knepley const PetscReal rp = (x[2]+1) * 0.5*(ro-ri) + ri; 1408b438e21SMatthew G. Knepley const PetscReal phip = PetscAtan2Real(x[1], x[0]); 1418b438e21SMatthew G. Knepley const PetscReal thetap = 0.5*PETSC_PI * (1.0 - ((((phip <= pi4) && (phip >= -pi4)) || ((phip >= 3.0*pi4) || (phip <= -3.0*pi4))) ? PetscAbsReal(x[0]) : PetscAbsReal(x[1]))); 1428b438e21SMatthew G. Knepley 1438b438e21SMatthew G. Knepley xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip); 1448b438e21SMatthew G. Knepley xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip); 1458b438e21SMatthew G. Knepley xp[2] = rp * PetscSinReal(thetap); 1468b438e21SMatthew G. Knepley } 1478b438e21SMatthew G. Knepley 14839290ff6SMatthew G. Knepley static PetscErrorCode DMCreateCoordinateDisc(DM dm) 1498b438e21SMatthew G. Knepley { 1508b438e21SMatthew G. Knepley DM cdm; 1518b438e21SMatthew G. Knepley PetscFE fe; 15239290ff6SMatthew G. Knepley DMPolytopeType ct; 15339290ff6SMatthew G. Knepley PetscInt dim, dE, cStart; 15439290ff6SMatthew G. Knepley PetscBool simplex; 1558b438e21SMatthew G. Knepley PetscErrorCode ierr; 1568b438e21SMatthew G. Knepley 1578b438e21SMatthew G. Knepley PetscFunctionBegin; 1588b438e21SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 1598b438e21SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1608b438e21SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr); 16139290ff6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, NULL);CHKERRQ(ierr); 16239290ff6SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 16339290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 16430602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe);CHKERRQ(ierr); 16539290ff6SMatthew G. Knepley ierr = DMProjectCoordinates(dm, fe);CHKERRQ(ierr); 1668b438e21SMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 1678b438e21SMatthew G. Knepley PetscFunctionReturn(0); 1688b438e21SMatthew G. Knepley } 1698b438e21SMatthew G. Knepley 1708b438e21SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 1718b438e21SMatthew G. Knepley { 1728b438e21SMatthew G. Knepley DM cdm; 1738b438e21SMatthew G. Knepley PetscDS cds; 1748b438e21SMatthew G. Knepley PetscErrorCode ierr; 1758b438e21SMatthew G. Knepley 1768b438e21SMatthew G. Knepley PetscFunctionBegin; 17730602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 17830602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1798b438e21SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 18030602db0SMatthew G. Knepley 18130602db0SMatthew G. Knepley if (ctx->coordSpace) {ierr = DMCreateCoordinateDisc(*dm);CHKERRQ(ierr);} 1828b438e21SMatthew G. Knepley switch (ctx->meshTransform) { 1838b438e21SMatthew G. Knepley case TRANSFORM_NONE: 1848b438e21SMatthew G. Knepley ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr); 1858b438e21SMatthew G. Knepley break; 1868b438e21SMatthew G. Knepley case TRANSFORM_SHEAR: 1878b438e21SMatthew G. Knepley ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr); 1888b438e21SMatthew G. Knepley break; 18960c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 19060c1a66aSMatthew G. Knepley ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 19160c1a66aSMatthew G. Knepley ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr); 19260c1a66aSMatthew G. Knepley ierr = PetscDSSetConstants(cds, 4, ctx->transformData);CHKERRQ(ierr); 19360c1a66aSMatthew G. Knepley ierr = DMPlexRemapGeometry(*dm, 0.0, f0_flare);CHKERRQ(ierr); 19460c1a66aSMatthew G. Knepley break; 1958b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 1968b438e21SMatthew G. Knepley ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 1978b438e21SMatthew G. Knepley ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr); 1988b438e21SMatthew G. Knepley ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr); 1998b438e21SMatthew G. Knepley ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr); 2008b438e21SMatthew G. Knepley break; 2018b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 2028b438e21SMatthew G. Knepley ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 2038b438e21SMatthew G. Knepley ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr); 2048b438e21SMatthew G. Knepley ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr); 2058b438e21SMatthew G. Knepley ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr); 2068b438e21SMatthew G. Knepley break; 207*98921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform); 2088b438e21SMatthew G. Knepley } 2098b438e21SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 2108b438e21SMatthew G. Knepley PetscFunctionReturn(0); 2118b438e21SMatthew G. Knepley } 2128b438e21SMatthew G. Knepley 2138b438e21SMatthew G. Knepley static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2148b438e21SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2158b438e21SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2168b438e21SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[]) 2178b438e21SMatthew G. Knepley { 2188b438e21SMatthew G. Knepley vol[0] = 1.; 2198b438e21SMatthew G. Knepley } 2208b438e21SMatthew G. Knepley 2218b438e21SMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx) 2228b438e21SMatthew G. Knepley { 2238b438e21SMatthew G. Knepley PetscDS ds; 2248b438e21SMatthew G. Knepley PetscFE fe; 22539290ff6SMatthew G. Knepley DMPolytopeType ct; 22639290ff6SMatthew G. Knepley PetscInt dim, cStart; 22739290ff6SMatthew G. Knepley PetscBool simplex; 2288b438e21SMatthew G. Knepley PetscErrorCode ierr; 2298b438e21SMatthew G. Knepley 2308b438e21SMatthew G. Knepley PetscFunctionBeginUser; 23139290ff6SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 23239290ff6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 23339290ff6SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 23439290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 23539290ff6SMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 2368b438e21SMatthew G. Knepley ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr); 2371e1ea65dSPierre Jolivet ierr = DMAddField(dm, NULL, (PetscObject) fe);CHKERRQ(ierr); 2388b438e21SMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 2398b438e21SMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 2408b438e21SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 2418b438e21SMatthew G. Knepley ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr); 2428b438e21SMatthew G. Knepley PetscFunctionReturn(0); 2438b438e21SMatthew G. Knepley } 2448b438e21SMatthew G. Knepley 2458b438e21SMatthew G. Knepley static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx) 2468b438e21SMatthew G. Knepley { 2478b438e21SMatthew G. Knepley Vec u; 2488b438e21SMatthew G. Knepley PetscScalar result; 2498b438e21SMatthew G. Knepley PetscReal vol, tol = ctx->tol; 2508b438e21SMatthew G. Knepley PetscErrorCode ierr; 2518b438e21SMatthew G. Knepley 2528b438e21SMatthew G. Knepley PetscFunctionBeginUser; 2538b438e21SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 2548b438e21SMatthew G. Knepley ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr); 2558b438e21SMatthew G. Knepley vol = PetscRealPart(result); 2568b438e21SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 2578b438e21SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr); 2588b438e21SMatthew G. Knepley if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) { 259*98921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Calculated volume %g != %g actual volume (error %g > %g tol)", (double) vol, (double) ctx->volume, (double) PetscAbsReal(ctx->volume - vol), (double) tol); 2608b438e21SMatthew G. Knepley } 2618b438e21SMatthew G. Knepley PetscFunctionReturn(0); 2628b438e21SMatthew G. Knepley } 2638b438e21SMatthew G. Knepley 2648b438e21SMatthew G. Knepley int main(int argc, char **argv) 2658b438e21SMatthew G. Knepley { 26639290ff6SMatthew G. Knepley DM dm; 2678b438e21SMatthew G. Knepley AppCtx user; 2688b438e21SMatthew G. Knepley PetscErrorCode ierr; 2698b438e21SMatthew G. Knepley 2708b438e21SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 2718b438e21SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 27239290ff6SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 27339290ff6SMatthew G. Knepley ierr = CreateDiscretization(dm, &user);CHKERRQ(ierr); 27439290ff6SMatthew G. Knepley ierr = CheckVolume(dm, &user);CHKERRQ(ierr); 27539290ff6SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 2768b438e21SMatthew G. Knepley ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr); 2778b438e21SMatthew G. Knepley ierr = PetscFree(user.transformData);CHKERRQ(ierr); 2788b438e21SMatthew G. Knepley ierr = PetscFinalize(); 2798b438e21SMatthew G. Knepley return ierr; 2808b438e21SMatthew G. Knepley } 2818b438e21SMatthew G. Knepley 2828b438e21SMatthew G. Knepley /*TEST 2838b438e21SMatthew G. Knepley 28430602db0SMatthew G. Knepley testset: 28530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0 28630602db0SMatthew G. Knepley 2878b438e21SMatthew G. Knepley test: 2888b438e21SMatthew G. Knepley suffix: square_0 28930602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 2908b438e21SMatthew G. Knepley 2918b438e21SMatthew G. Knepley test: 2928b438e21SMatthew G. Knepley suffix: square_1 29330602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 2948b438e21SMatthew G. Knepley 2958b438e21SMatthew G. Knepley test: 2968b438e21SMatthew G. Knepley suffix: square_2 29730602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 2988b438e21SMatthew G. Knepley 2998b438e21SMatthew G. Knepley test: 3008b438e21SMatthew G. Knepley suffix: square_3 30130602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 30230602db0SMatthew G. Knepley 30330602db0SMatthew G. Knepley testset: 30430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0 3058b438e21SMatthew G. Knepley 3068b438e21SMatthew G. Knepley test: 3078b438e21SMatthew G. Knepley suffix: cube_0 30830602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 3098b438e21SMatthew G. Knepley 3108b438e21SMatthew G. Knepley test: 3118b438e21SMatthew G. Knepley suffix: cube_1 31230602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 3138b438e21SMatthew G. Knepley 3148b438e21SMatthew G. Knepley test: 3158b438e21SMatthew G. Knepley suffix: cube_2 31630602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 3178b438e21SMatthew G. Knepley 3188b438e21SMatthew G. Knepley test: 3198b438e21SMatthew G. Knepley suffix: cube_3 32030602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 32130602db0SMatthew G. Knepley 32230602db0SMatthew G. Knepley testset: 32330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0 3248b438e21SMatthew G. Knepley 3258b438e21SMatthew G. Knepley test: 3268b438e21SMatthew G. Knepley suffix: shear_0 32730602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3288b438e21SMatthew G. Knepley 3298b438e21SMatthew G. Knepley test: 3308b438e21SMatthew G. Knepley suffix: shear_1 33130602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3328b438e21SMatthew G. Knepley 3338b438e21SMatthew G. Knepley test: 3348b438e21SMatthew G. Knepley suffix: shear_2 33530602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3368b438e21SMatthew G. Knepley 3378b438e21SMatthew G. Knepley test: 3388b438e21SMatthew G. Knepley suffix: shear_3 33930602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 34030602db0SMatthew G. Knepley 34130602db0SMatthew G. Knepley testset: 34230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0 3438b438e21SMatthew G. Knepley 3448b438e21SMatthew G. Knepley test: 3458b438e21SMatthew G. Knepley suffix: shear_4 34630602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3478b438e21SMatthew G. Knepley 3488b438e21SMatthew G. Knepley test: 3498b438e21SMatthew G. Knepley suffix: shear_5 35030602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3518b438e21SMatthew G. Knepley 3528b438e21SMatthew G. Knepley test: 3538b438e21SMatthew G. Knepley suffix: shear_6 35430602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 3558b438e21SMatthew G. Knepley 3568b438e21SMatthew G. Knepley test: 3578b438e21SMatthew G. Knepley suffix: shear_7 35830602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 35930602db0SMatthew G. Knepley 36030602db0SMatthew G. Knepley testset: 36160c1a66aSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower 1.,-1. -dm_plex_box_upper 3.,1. \ 36260c1a66aSMatthew G. Knepley -dm_coord_petscspace_degree 1 -mesh_transform flare -volume 8. 36360c1a66aSMatthew G. Knepley 36460c1a66aSMatthew G. Knepley test: 36560c1a66aSMatthew G. Knepley suffix: flare_0 36660c1a66aSMatthew G. Knepley args: 36760c1a66aSMatthew G. Knepley 36860c1a66aSMatthew G. Knepley test: 36960c1a66aSMatthew G. Knepley suffix: flare_1 37060c1a66aSMatthew G. Knepley args: -dm_refine 2 37160c1a66aSMatthew G. Knepley 37260c1a66aSMatthew G. Knepley testset: 37330602db0SMatthew G. Knepley # Area: 3/4 \pi = 2.3562 37430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0 3758b438e21SMatthew G. Knepley 3768b438e21SMatthew G. Knepley test: 3778b438e21SMatthew G. Knepley # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2 3788b438e21SMatthew G. Knepley suffix: annulus_0 3798b438e21SMatthew G. Knepley requires: double 38030602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -volume 1.5 3818b438e21SMatthew G. Knepley 3828b438e21SMatthew G. Knepley test: 3838b438e21SMatthew G. Knepley suffix: annulus_1 3848b438e21SMatthew G. Knepley requires: double 38530602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016 3868b438e21SMatthew G. Knepley 3878b438e21SMatthew G. Knepley test: 3888b438e21SMatthew G. Knepley suffix: annulus_2 3898b438e21SMatthew G. Knepley requires: double 39030602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038 3918b438e21SMatthew G. Knepley 3928b438e21SMatthew G. Knepley test: 3938b438e21SMatthew G. Knepley suffix: annulus_3 3948b438e21SMatthew G. Knepley requires: double 39530602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6 3968b438e21SMatthew G. Knepley 3978b438e21SMatthew G. Knepley test: 3988b438e21SMatthew G. Knepley suffix: annulus_4 3998b438e21SMatthew G. Knepley requires: double 40030602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012 4018b438e21SMatthew G. Knepley 4028b438e21SMatthew G. Knepley test: 4038b438e21SMatthew G. Knepley suffix: annulus_5 4048b438e21SMatthew G. Knepley requires: double 40530602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7 40630602db0SMatthew G. Knepley 40730602db0SMatthew G. Knepley testset: 40830602db0SMatthew G. Knepley # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238 40930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -mesh_transform shell -volume 14.66076571675238 -dm_coord_space 0 4108b438e21SMatthew G. Knepley 4118b438e21SMatthew G. Knepley test: 4128b438e21SMatthew G. Knepley suffix: shell_0 4138b438e21SMatthew G. Knepley requires: double 41430602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7 4158b438e21SMatthew G. Knepley 4168b438e21SMatthew G. Knepley test: 4178b438e21SMatthew G. Knepley suffix: shell_1 4188b438e21SMatthew G. Knepley requires: double 41930602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1 4208b438e21SMatthew G. Knepley 4218b438e21SMatthew G. Knepley test: 4228b438e21SMatthew G. Knepley suffix: shell_2 4238b438e21SMatthew G. Knepley requires: double 42430602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1 4258b438e21SMatthew G. Knepley 4268b438e21SMatthew G. Knepley test: 4278b438e21SMatthew G. Knepley suffix: shell_3 4288b438e21SMatthew G. Knepley requires: double 42930602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02 4308b438e21SMatthew G. Knepley 4318b438e21SMatthew G. Knepley test: 4328b438e21SMatthew G. Knepley suffix: shell_4 4338b438e21SMatthew G. Knepley requires: double 43430602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006 43539290ff6SMatthew G. Knepley 43639290ff6SMatthew G. Knepley test: 43739290ff6SMatthew G. Knepley # Volume: 1.0 43839290ff6SMatthew G. Knepley suffix: gmsh_q2 43939290ff6SMatthew G. Knepley requires: double 44030602db0SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 44139290ff6SMatthew G. Knepley 44239290ff6SMatthew G. Knepley test: 44339290ff6SMatthew G. Knepley # Volume: 1.0 44439290ff6SMatthew G. Knepley suffix: gmsh_q3 44539290ff6SMatthew G. Knepley requires: double 44630602db0SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q3.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 4478b438e21SMatthew G. Knepley 4488b438e21SMatthew G. Knepley TEST*/ 449