xref: /petsc/src/dm/impls/plex/tutorials/ex9.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
10f108dbd7SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11c4762a1bSJed Brown {
12f108dbd7SJacob Faibussowitsch   PetscFunctionBeginUser;
13c4762a1bSJed Brown   options->report      = PETSC_FALSE;
14f108dbd7SJacob Faibussowitsch   options->tol         = 0.5;
15c4762a1bSJed Brown   options->condLimit   = PETSC_DETERMINE;
16c4762a1bSJed Brown 
17*d0609cedSBarry 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));
21*d0609cedSBarry Smith   PetscOptionsEnd();
22f108dbd7SJacob Faibussowitsch   PetscFunctionReturn(0);
23f108dbd7SJacob Faibussowitsch }
24c4762a1bSJed Brown 
25f108dbd7SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
26f108dbd7SJacob 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"));
32c4762a1bSJed Brown   PetscFunctionReturn(0);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35c4762a1bSJed Brown int main(int argc, char **argv)
36c4762a1bSJed Brown {
37c4762a1bSJed Brown   DM             dm;
38f108dbd7SJacob Faibussowitsch   DMLabel        OQLabel;
39f108dbd7SJacob Faibussowitsch   Vec            OQ;
40c4762a1bSJed Brown   AppCtx         ctx;
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
439566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
449566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
459566063dSJacob Faibussowitsch   PetscCall(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel));
479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&OQ));
489566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
50b122ec5aSJacob Faibussowitsch   return 0;
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /*TEST
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   test:
56c4762a1bSJed Brown     suffix: 0
57c4762a1bSJed Brown     requires: exodusii
58c4762a1bSJed Brown     nsize: {{1 2}}
59e600fa54SMatthew G. Knepley     args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   test:
62c4762a1bSJed Brown     suffix: 1
6330602db0SMatthew G. Knepley     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report
64c4762a1bSJed Brown 
65f108dbd7SJacob Faibussowitsch   testset:
66f108dbd7SJacob Faibussowitsch     args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view
67f108dbd7SJacob Faibussowitsch 
68f108dbd7SJacob Faibussowitsch     test:
69f108dbd7SJacob Faibussowitsch       suffix: box_1
70f108dbd7SJacob Faibussowitsch       nsize: 1
7130602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
72f108dbd7SJacob Faibussowitsch 
73f108dbd7SJacob Faibussowitsch     test:
74f108dbd7SJacob Faibussowitsch       suffix: box_2
75f108dbd7SJacob Faibussowitsch       nsize: 2
76e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0
77f108dbd7SJacob Faibussowitsch 
78f108dbd7SJacob Faibussowitsch     test:
79f108dbd7SJacob Faibussowitsch       suffix: mesh_1
80f108dbd7SJacob Faibussowitsch       nsize: 1
816ed19f2fSJacob Faibussowitsch       requires: exodusii
8230602db0SMatthew G. Knepley       args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
83f108dbd7SJacob Faibussowitsch 
84f108dbd7SJacob Faibussowitsch     test:
85f108dbd7SJacob Faibussowitsch       suffix: mesh_2
86f108dbd7SJacob Faibussowitsch       nsize: 2
876ed19f2fSJacob Faibussowitsch       requires: exodusii
88e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
89c4762a1bSJed Brown TEST*/
90