1 static char help[] = "Tests for high order geometry\n\n"; 2 3 #include <petscdmplex.h> 4 #include <petscds.h> 5 6 typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform; 7 const char * const TransformTypes[] = {"none", "shear", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL}; 8 9 typedef struct { 10 PetscBool coordSpace; /* Flag to create coordinate space */ 11 Transform meshTransform; /* Transform for initial box mesh */ 12 PetscReal *transformDataReal; /* Parameters for mesh transform */ 13 PetscScalar *transformData; /* Parameters for mesh transform */ 14 PetscReal volume; /* Analytical volume of the mesh */ 15 PetscReal tol; /* Tolerance for volume check */ 16 } AppCtx; 17 18 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 19 { 20 PetscInt n = 0, i; 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 options->coordSpace = PETSC_TRUE; 25 options->meshTransform = TRANSFORM_NONE; 26 options->transformDataReal = NULL; 27 options->transformData = NULL; 28 options->volume = -1.0; 29 options->tol = PETSC_SMALL; 30 31 ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr); 32 ierr = PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL);CHKERRQ(ierr); 33 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); 34 switch (options->meshTransform) { 35 case TRANSFORM_NONE: break; 36 case TRANSFORM_SHEAR: 37 n = 2; 38 ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr); 39 for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0; 40 ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr); 41 break; 42 case TRANSFORM_ANNULUS: 43 n = 2; 44 ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr); 45 options->transformData[0] = 1.0; 46 options->transformData[1] = 2.0; 47 ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr); 48 break; 49 case TRANSFORM_SHELL: 50 n = 2; 51 ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr); 52 options->transformData[0] = 1.0; 53 options->transformData[1] = 2.0; 54 ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr); 55 break; 56 default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform); 57 } 58 ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr); 59 ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr); 60 ierr = PetscOptionsEnd(); 61 PetscFunctionReturn(0); 62 } 63 64 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 65 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 66 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 67 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 68 { 69 const PetscInt Nc = uOff[1] - uOff[0]; 70 PetscInt c; 71 72 for (c = 0; c < Nc; ++c) f0[c] = u[c]; 73 } 74 75 /* 76 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 77 will correspond to the top and bottom of our square. So 78 79 (0,0)--(1,0) ==> (1,0)--(2,0) Just a shift of (1,0) 80 (0,1)--(1,1) ==> (0,1)--(0,2) Switch x and y 81 82 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: 83 84 (x, y) ==> (x+1, \pi/2 y) in (r', \theta') space 85 ==> ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space 86 */ 87 static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, 88 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 89 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 90 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 91 { 92 const PetscReal ri = PetscRealPart(constants[0]); 93 const PetscReal ro = PetscRealPart(constants[1]); 94 95 xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]); 96 xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]); 97 } 98 99 /* 100 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 101 lower hemisphere and the upper surface onto the top, letting z be the radius. 102 103 (x, y) ==> ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x) in (r', \theta', \phi') space 104 ==> ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space 105 */ 106 static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux, 107 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 108 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 109 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 110 { 111 const PetscReal pi4 = PETSC_PI/4.0; 112 const PetscReal ri = PetscRealPart(constants[0]); 113 const PetscReal ro = PetscRealPart(constants[1]); 114 const PetscReal rp = (x[2]+1) * 0.5*(ro-ri) + ri; 115 const PetscReal phip = PetscAtan2Real(x[1], x[0]); 116 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]))); 117 118 xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip); 119 xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip); 120 xp[2] = rp * PetscSinReal(thetap); 121 } 122 123 static PetscErrorCode DMCreateCoordinateDisc(DM dm) 124 { 125 DM cdm; 126 PetscFE fe; 127 DMPolytopeType ct; 128 PetscInt dim, dE, cStart; 129 PetscBool simplex; 130 PetscErrorCode ierr; 131 132 PetscFunctionBegin; 133 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 134 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 135 ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr); 136 ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, NULL);CHKERRQ(ierr); 137 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 138 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 139 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe);CHKERRQ(ierr); 140 ierr = DMProjectCoordinates(dm, fe);CHKERRQ(ierr); 141 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 142 PetscFunctionReturn(0); 143 } 144 145 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 146 { 147 DM cdm; 148 PetscDS cds; 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 153 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 154 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 155 156 if (ctx->coordSpace) {ierr = DMCreateCoordinateDisc(*dm);CHKERRQ(ierr);} 157 switch (ctx->meshTransform) { 158 case TRANSFORM_NONE: 159 ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr); 160 break; 161 case TRANSFORM_SHEAR: 162 ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr); 163 break; 164 case TRANSFORM_ANNULUS: 165 ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 166 ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr); 167 ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr); 168 ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr); 169 break; 170 case TRANSFORM_SHELL: 171 ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 172 ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr); 173 ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr); 174 ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr); 175 break; 176 default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform); 177 } 178 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 179 PetscFunctionReturn(0); 180 } 181 182 static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux, 183 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 184 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 185 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[]) 186 { 187 vol[0] = 1.; 188 } 189 190 static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx) 191 { 192 PetscDS ds; 193 PetscFE fe; 194 DMPolytopeType ct; 195 PetscInt dim, cStart; 196 PetscBool simplex; 197 PetscErrorCode ierr; 198 199 PetscFunctionBeginUser; 200 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 201 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 202 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 203 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 204 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 205 ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr); 206 ierr = DMAddField(dm, NULL, (PetscObject) fe); 207 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 208 ierr = DMCreateDS(dm);CHKERRQ(ierr); 209 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 210 ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr); 211 PetscFunctionReturn(0); 212 } 213 214 static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx) 215 { 216 Vec u; 217 PetscScalar result; 218 PetscReal vol, tol = ctx->tol; 219 PetscErrorCode ierr; 220 221 PetscFunctionBeginUser; 222 ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 223 ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr); 224 vol = PetscRealPart(result); 225 ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 226 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr); 227 if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) { 228 SETERRQ4(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); 229 } 230 PetscFunctionReturn(0); 231 } 232 233 int main(int argc, char **argv) 234 { 235 DM dm; 236 AppCtx user; 237 PetscErrorCode ierr; 238 239 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 240 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 241 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 242 ierr = CreateDiscretization(dm, &user);CHKERRQ(ierr); 243 ierr = CheckVolume(dm, &user);CHKERRQ(ierr); 244 ierr = DMDestroy(&dm);CHKERRQ(ierr); 245 ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr); 246 ierr = PetscFree(user.transformData);CHKERRQ(ierr); 247 ierr = PetscFinalize(); 248 return ierr; 249 } 250 251 /*TEST 252 253 testset: 254 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 255 256 test: 257 suffix: square_0 258 args: -dm_coord_petscspace_degree 1 259 260 test: 261 suffix: square_1 262 args: -dm_coord_petscspace_degree 2 263 264 test: 265 suffix: square_2 266 args: -dm_refine 1 -dm_coord_petscspace_degree 1 267 268 test: 269 suffix: square_3 270 args: -dm_refine 1 -dm_coord_petscspace_degree 2 271 272 testset: 273 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 274 275 test: 276 suffix: cube_0 277 args: -dm_coord_petscspace_degree 1 278 279 test: 280 suffix: cube_1 281 args: -dm_coord_petscspace_degree 2 282 283 test: 284 suffix: cube_2 285 args: -dm_refine 1 -dm_coord_petscspace_degree 1 286 287 test: 288 suffix: cube_3 289 args: -dm_refine 1 -dm_coord_petscspace_degree 2 290 291 testset: 292 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 293 294 test: 295 suffix: shear_0 296 args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 297 298 test: 299 suffix: shear_1 300 args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 301 302 test: 303 suffix: shear_2 304 args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 305 306 test: 307 suffix: shear_3 308 args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 309 310 testset: 311 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 312 313 test: 314 suffix: shear_4 315 args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 316 317 test: 318 suffix: shear_5 319 args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 320 321 test: 322 suffix: shear_6 323 args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 324 325 test: 326 suffix: shear_7 327 args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 328 329 testset: 330 # Area: 3/4 \pi = 2.3562 331 args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0 332 333 test: 334 # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2 335 suffix: annulus_0 336 requires: double 337 args: -dm_coord_petscspace_degree 1 -volume 1.5 338 339 test: 340 suffix: annulus_1 341 requires: double 342 args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016 343 344 test: 345 suffix: annulus_2 346 requires: double 347 args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038 348 349 test: 350 suffix: annulus_3 351 requires: double 352 args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6 353 354 test: 355 suffix: annulus_4 356 requires: double 357 args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012 358 359 test: 360 suffix: annulus_5 361 requires: double 362 args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7 363 364 testset: 365 # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238 366 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 367 368 test: 369 suffix: shell_0 370 requires: double 371 args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7 372 373 test: 374 suffix: shell_1 375 requires: double 376 args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1 377 378 test: 379 suffix: shell_2 380 requires: double 381 args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1 382 383 test: 384 suffix: shell_3 385 requires: double 386 args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02 387 388 test: 389 suffix: shell_4 390 requires: double 391 args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006 392 393 test: 394 # Volume: 1.0 395 suffix: gmsh_q2 396 requires: double 397 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 398 399 test: 400 # Volume: 1.0 401 suffix: gmsh_q3 402 requires: double 403 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 404 405 TEST*/ 406