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