1c4762a1bSJed Brown static char help[] = "Evaluate the shape quality of a mesh\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown PetscBool report; /* Print a quality report */
7f108dbd7SJacob Faibussowitsch PetscReal condLimit, tol; /* Condition number limit for cell output */
8c4762a1bSJed Brown } AppCtx;
9c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)10d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11d71ae5a4SJacob Faibussowitsch {
12f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser;
13c4762a1bSJed Brown options->report = PETSC_FALSE;
14f108dbd7SJacob Faibussowitsch options->tol = 0.5;
15c4762a1bSJed Brown options->condLimit = PETSC_DETERMINE;
16c4762a1bSJed Brown
17d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX");
189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL));
199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL));
209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL));
21d0609cedSBarry Smith PetscOptionsEnd();
223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
23f108dbd7SJacob Faibussowitsch }
24c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)25d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
26d71ae5a4SJacob Faibussowitsch {
27f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser;
289566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
299566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
319566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown
main(int argc,char ** argv)35d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
36d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown DM dm;
38f108dbd7SJacob Faibussowitsch DMLabel OQLabel;
39f108dbd7SJacob Faibussowitsch Vec OQ;
40c4762a1bSJed Brown AppCtx ctx;
41c4762a1bSJed Brown
42327415f7SBarry Smith PetscFunctionBeginUser;
439566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
449566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
459566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
468fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm));
479566063dSJacob Faibussowitsch PetscCall(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit));
489566063dSJacob Faibussowitsch PetscCall(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel));
499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&OQ));
509566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
519566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
52b122ec5aSJacob Faibussowitsch return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown
57c4762a1bSJed Brown test:
58c4762a1bSJed Brown suffix: 0
59c4762a1bSJed Brown requires: exodusii
60c4762a1bSJed Brown nsize: {{1 2}}
61*46ac1a18SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report
62c4762a1bSJed Brown
63c4762a1bSJed Brown test:
64c4762a1bSJed Brown suffix: 1
6530602db0SMatthew G. Knepley args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report
66c4762a1bSJed Brown
67f108dbd7SJacob Faibussowitsch testset:
68f108dbd7SJacob Faibussowitsch args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view
69f108dbd7SJacob Faibussowitsch
70f108dbd7SJacob Faibussowitsch test:
71f108dbd7SJacob Faibussowitsch suffix: box_1
72f108dbd7SJacob Faibussowitsch nsize: 1
7330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
74f108dbd7SJacob Faibussowitsch
75f108dbd7SJacob Faibussowitsch test:
76f108dbd7SJacob Faibussowitsch suffix: box_2
77f108dbd7SJacob Faibussowitsch nsize: 2
78e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
79f108dbd7SJacob Faibussowitsch
80f108dbd7SJacob Faibussowitsch test:
81f108dbd7SJacob Faibussowitsch suffix: mesh_1
82f108dbd7SJacob Faibussowitsch nsize: 1
836ed19f2fSJacob Faibussowitsch requires: exodusii
8430602db0SMatthew G. Knepley args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
85f108dbd7SJacob Faibussowitsch
86f108dbd7SJacob Faibussowitsch test:
87f108dbd7SJacob Faibussowitsch suffix: mesh_2
88f108dbd7SJacob Faibussowitsch nsize: 2
896ed19f2fSJacob Faibussowitsch requires: exodusii
90e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
91c4762a1bSJed Brown TEST*/
92