1 static char help[] = "Evaluate the shape quality of a mesh\n\n"; 2 3 #include <petscdmplex.h> 4 5 typedef struct { 6 PetscBool report; /* Print a quality report */ 7 PetscReal condLimit, tol; /* Condition number limit for cell output */ 8 } AppCtx; 9 10 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11 { 12 PetscFunctionBeginUser; 13 options->report = PETSC_FALSE; 14 options->tol = 0.5; 15 options->condLimit = PETSC_DETERMINE; 16 17 PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX"); 18 PetscCall(PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL)); 19 PetscCall(PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL)); 20 PetscCall(PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL)); 21 PetscOptionsEnd(); 22 PetscFunctionReturn(0); 23 } 24 25 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 26 { 27 PetscFunctionBeginUser; 28 PetscCall(DMCreate(comm, dm)); 29 PetscCall(DMSetType(*dm, DMPLEX)); 30 PetscCall(DMSetFromOptions(*dm)); 31 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 32 PetscFunctionReturn(0); 33 } 34 35 int main(int argc, char **argv) 36 { 37 DM dm; 38 DMLabel OQLabel; 39 Vec OQ; 40 AppCtx ctx; 41 42 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 43 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 44 PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 45 PetscCall(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit)); 46 PetscCall(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel)); 47 PetscCall(VecDestroy(&OQ)); 48 PetscCall(DMDestroy(&dm)); 49 PetscCall(PetscFinalize()); 50 return 0; 51 } 52 53 /*TEST 54 55 test: 56 suffix: 0 57 requires: exodusii 58 nsize: {{1 2}} 59 args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report 60 61 test: 62 suffix: 1 63 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report 64 65 testset: 66 args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view 67 68 test: 69 suffix: box_1 70 nsize: 1 71 args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 72 73 test: 74 suffix: box_2 75 nsize: 2 76 args: -petscpartitioner_type simple -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 77 78 test: 79 suffix: mesh_1 80 nsize: 1 81 requires: exodusii 82 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 83 84 test: 85 suffix: mesh_2 86 nsize: 2 87 requires: exodusii 88 args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 89 TEST*/ 90