1c4762a1bSJed Brown static char help[] = "Evaluate the shape quality of a mesh\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown PetscBool report; /* Print a quality report */ 7f108dbd7SJacob Faibussowitsch PetscReal condLimit, tol; /* Condition number limit for cell output */ 8c4762a1bSJed Brown } AppCtx; 9c4762a1bSJed Brown 10f108dbd7SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11c4762a1bSJed Brown { 12f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser; 13c4762a1bSJed Brown options->report = PETSC_FALSE; 14f108dbd7SJacob Faibussowitsch options->tol = 0.5; 15c4762a1bSJed Brown options->condLimit = PETSC_DETERMINE; 16c4762a1bSJed Brown 17d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX"); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL)); 21d0609cedSBarry Smith PetscOptionsEnd(); 22f108dbd7SJacob Faibussowitsch PetscFunctionReturn(0); 23f108dbd7SJacob Faibussowitsch } 24c4762a1bSJed Brown 25f108dbd7SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 26f108dbd7SJacob Faibussowitsch { 27f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser; 289566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 299566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 319566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35c4762a1bSJed Brown int main(int argc, char **argv) 36c4762a1bSJed Brown { 37c4762a1bSJed Brown DM dm; 38f108dbd7SJacob Faibussowitsch DMLabel OQLabel; 39f108dbd7SJacob Faibussowitsch Vec OQ; 40c4762a1bSJed Brown AppCtx ctx; 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 439566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 449566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 45*8fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 469566063dSJacob Faibussowitsch PetscCall(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit)); 479566063dSJacob Faibussowitsch PetscCall(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel)); 489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&OQ)); 499566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 51b122ec5aSJacob Faibussowitsch return 0; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /*TEST 55c4762a1bSJed Brown 56c4762a1bSJed Brown test: 57c4762a1bSJed Brown suffix: 0 58c4762a1bSJed Brown requires: exodusii 59c4762a1bSJed Brown nsize: {{1 2}} 60e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report 61c4762a1bSJed Brown 62c4762a1bSJed Brown test: 63c4762a1bSJed Brown suffix: 1 6430602db0SMatthew G. Knepley args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report 65c4762a1bSJed Brown 66f108dbd7SJacob Faibussowitsch testset: 67f108dbd7SJacob Faibussowitsch args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view 68f108dbd7SJacob Faibussowitsch 69f108dbd7SJacob Faibussowitsch test: 70f108dbd7SJacob Faibussowitsch suffix: box_1 71f108dbd7SJacob Faibussowitsch nsize: 1 7230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 73f108dbd7SJacob Faibussowitsch 74f108dbd7SJacob Faibussowitsch test: 75f108dbd7SJacob Faibussowitsch suffix: box_2 76f108dbd7SJacob Faibussowitsch nsize: 2 77e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 78f108dbd7SJacob Faibussowitsch 79f108dbd7SJacob Faibussowitsch test: 80f108dbd7SJacob Faibussowitsch suffix: mesh_1 81f108dbd7SJacob Faibussowitsch nsize: 1 826ed19f2fSJacob Faibussowitsch requires: exodusii 8330602db0SMatthew G. Knepley args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 84f108dbd7SJacob Faibussowitsch 85f108dbd7SJacob Faibussowitsch test: 86f108dbd7SJacob Faibussowitsch suffix: mesh_2 87f108dbd7SJacob Faibussowitsch nsize: 2 886ed19f2fSJacob Faibussowitsch requires: exodusii 89e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 90c4762a1bSJed Brown TEST*/ 91