16e160c9fSMatthew G. Knepley static char help[] = "Tests for geometry models\n\n"; 26e160c9fSMatthew G. Knepley 36e160c9fSMatthew G. Knepley #include <petscdmplex.h> 46e160c9fSMatthew G. Knepley #include <petscds.h> 56e160c9fSMatthew G. Knepley 66e160c9fSMatthew G. Knepley typedef struct { 76e160c9fSMatthew G. Knepley PetscReal exactVol; // The exact volume of the shape 86e160c9fSMatthew G. Knepley PetscReal volTol; // The relative tolerance for checking the volume 96e160c9fSMatthew G. Knepley } AppCtx; 106e160c9fSMatthew G. Knepley 116e160c9fSMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 126e160c9fSMatthew G. Knepley { 136e160c9fSMatthew G. Knepley f0[0] = 1.0; 146e160c9fSMatthew G. Knepley } 156e160c9fSMatthew G. Knepley 166e160c9fSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 176e160c9fSMatthew G. Knepley { 186e160c9fSMatthew G. Knepley PetscFunctionBegin; 196e160c9fSMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 206e160c9fSMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 216e160c9fSMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 226e160c9fSMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 236e160c9fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 246e160c9fSMatthew G. Knepley } 256e160c9fSMatthew G. Knepley 266e160c9fSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 276e160c9fSMatthew G. Knepley { 286e160c9fSMatthew G. Knepley options->exactVol = 4. * PETSC_PI / 3.; 296e160c9fSMatthew G. Knepley options->volTol = PETSC_SMALL; 306e160c9fSMatthew G. Knepley 316e160c9fSMatthew G. Knepley PetscFunctionBeginUser; 326e160c9fSMatthew G. Knepley PetscOptionsBegin(comm, "", "Geometry Model Test Options", "DMPLEX"); 336e160c9fSMatthew G. Knepley PetscCall(PetscOptionsReal("-exact_vol", "Exact volume of the shape", __FILE__, options->exactVol, &options->exactVol, NULL)); 346e160c9fSMatthew G. Knepley PetscCall(PetscOptionsReal("-vol_tol", "Relative tolerance for checking the volume", __FILE__, options->volTol, &options->volTol, NULL)); 356e160c9fSMatthew G. Knepley PetscOptionsEnd(); 366e160c9fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 376e160c9fSMatthew G. Knepley } 386e160c9fSMatthew G. Knepley 396e160c9fSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm) 406e160c9fSMatthew G. Knepley { 416e160c9fSMatthew G. Knepley PetscFE fe; 426e160c9fSMatthew G. Knepley PetscDS ds; 436e160c9fSMatthew G. Knepley DMPolytopeType ct; 446e160c9fSMatthew G. Knepley PetscInt dim, cStart; 456e160c9fSMatthew G. Knepley 466e160c9fSMatthew G. Knepley PetscFunctionBeginUser; 476e160c9fSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 486e160c9fSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 496e160c9fSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 506e160c9fSMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, NULL, -1, &fe)); 516e160c9fSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 526e160c9fSMatthew G. Knepley PetscCall(DMCreateDS(dm)); 536e160c9fSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 546e160c9fSMatthew G. Knepley PetscCall(PetscDSSetObjective(ds, 0, identity)); 556e160c9fSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe)); 566e160c9fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 576e160c9fSMatthew G. Knepley } 586e160c9fSMatthew G. Knepley 596e160c9fSMatthew G. Knepley int main(int argc, char **argv) 606e160c9fSMatthew G. Knepley { 616e160c9fSMatthew G. Knepley DM dm; 626e160c9fSMatthew G. Knepley Vec u; 636e160c9fSMatthew G. Knepley PetscScalar volume; 646e160c9fSMatthew G. Knepley AppCtx user; 656e160c9fSMatthew G. Knepley 666e160c9fSMatthew G. Knepley PetscFunctionBeginUser; 676e160c9fSMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 686e160c9fSMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 696e160c9fSMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 706e160c9fSMatthew G. Knepley PetscCall(CreateDiscretization(dm)); 716e160c9fSMatthew G. Knepley PetscCall(DMGetGlobalVector(dm, &u)); 726e160c9fSMatthew G. Knepley PetscCall(VecSet(u, 0.)); 736e160c9fSMatthew G. Knepley PetscCall(DMPlexComputeIntegralFEM(dm, u, &volume, NULL)); 746e160c9fSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dm, &u)); 756e160c9fSMatthew G. Knepley PetscCheck(PetscAbsScalar((volume - user.exactVol) / user.exactVol) < user.volTol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid volume %g != %g", (double)PetscRealPart(volume), (double)user.exactVol); 766e160c9fSMatthew G. Knepley PetscCall(DMDestroy(&dm)); 776e160c9fSMatthew G. Knepley PetscCall(PetscFinalize()); 786e160c9fSMatthew G. Knepley return 0; 796e160c9fSMatthew G. Knepley } 806e160c9fSMatthew G. Knepley 816e160c9fSMatthew G. Knepley /*TEST 826e160c9fSMatthew G. Knepley 836e160c9fSMatthew G. Knepley # -dm_refine 6 -vol_tol 6e-4 works 846e160c9fSMatthew G. Knepley test: 856e160c9fSMatthew G. Knepley suffix: ball_0 866e160c9fSMatthew G. Knepley requires: ctetgen 876e160c9fSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_shape ball -dm_refine 3 -dm_geom_model ball -vol_tol 5e-2 88*3886731fSPierre Jolivet output_file: output/empty.out 896e160c9fSMatthew G. Knepley 906e160c9fSMatthew G. Knepley # -dm_refine 4 -vol_tol 2e-3 works 916e160c9fSMatthew G. Knepley test: 926e160c9fSMatthew G. Knepley suffix: cylinder_0 936e160c9fSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_shape cylinder -dm_plex_cylinder_num_refine 1 \ 946e160c9fSMatthew G. Knepley -dm_refine 1 -dm_geom_model cylinder -exact_vol 3.141592653589 -vol_tol 1e-1 95*3886731fSPierre Jolivet output_file: output/empty.out 966e160c9fSMatthew G. Knepley 976e160c9fSMatthew G. Knepley TEST*/ 98