Actual source code: ex9.c
1: static char help[] = "Evaluate the shape quality of a mesh\n\n";
3: #include <petscdmplex.h>
5: typedef struct {
6: PetscBool report; /* Print a quality report */
7: PetscReal condLimit, tol; /* Condition number limit for cell output */
8: } AppCtx;
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;
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(PETSC_SUCCESS);
23: }
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(PETSC_SUCCESS);
33: }
35: int main(int argc, char **argv)
36: {
37: DM dm;
38: DMLabel OQLabel;
39: Vec OQ;
40: AppCtx ctx;
42: PetscFunctionBeginUser;
43: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
44: PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
45: PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
46: PetscCall(DMGetCoordinatesLocalSetUp(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: }
55: /*TEST
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
63: test:
64: suffix: 1
65: args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report
67: testset:
68: args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view
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
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
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
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*/