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