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: {
15: options->report = PETSC_FALSE;
16: options->tol = 0.5;
17: options->condLimit = PETSC_DETERMINE;
19: PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX");
20: PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL);
21: PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL);
22: PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL);
23: PetscOptionsEnd();
24: return(0);
25: }
27: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
28: {
32: DMCreate(comm, dm);
33: DMSetType(*dm, DMPLEX);
34: DMSetFromOptions(*dm);
35: DMViewFromOptions(*dm, NULL, "-dm_view");
36: return(0);
37: }
39: int main(int argc, char **argv)
40: {
41: DM dm;
42: DMLabel OQLabel;
43: Vec OQ;
44: AppCtx ctx;
47: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
48: ProcessOptions(PETSC_COMM_WORLD, &ctx);
49: CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);
50: DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit);
51: DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel);
52: VecDestroy(&OQ);
53: DMDestroy(&dm);
54: PetscFinalize();
55: return ierr;
56: }
58: /*TEST
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
66: test:
67: suffix: 1
68: args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report
70: testset:
71: args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view
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
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
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
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*/