xref: /petsc/src/dm/impls/plex/tutorials/ex9.c (revision 99a7f59ec9f9ccee3cdf9fd88ce16d7f72105d85)
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(DMGetCoordinatesLocalSetUp(dm));
46   PetscCall(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit));
47   PetscCall(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel));
48   PetscCall(VecDestroy(&OQ));
49   PetscCall(DMDestroy(&dm));
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55 
56   test:
57     suffix: 0
58     requires: exodusii
59     nsize: {{1 2}}
60     args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report
61 
62   test:
63     suffix: 1
64     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report
65 
66   testset:
67     args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view
68 
69     test:
70       suffix: box_1
71       nsize: 1
72       args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
73 
74     test:
75       suffix: box_2
76       nsize: 2
77       args: -petscpartitioner_type simple -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
78 
79     test:
80       suffix: mesh_1
81       nsize: 1
82       requires: exodusii
83       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
84 
85     test:
86       suffix: mesh_2
87       nsize: 2
88       requires: exodusii
89       args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
90 TEST*/
91