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