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