xref: /petsc/src/dm/impls/plex/tutorials/ex9.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 {
12c4762a1bSJed Brown   PetscErrorCode ierr;
13c4762a1bSJed Brown 
14f108dbd7SJacob Faibussowitsch   PetscFunctionBeginUser;
15c4762a1bSJed Brown   options->report      = PETSC_FALSE;
16f108dbd7SJacob Faibussowitsch   options->tol         = 0.5;
17c4762a1bSJed Brown   options->condLimit   = PETSC_DETERMINE;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX");CHKERRQ(ierr);
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL));
23c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
24f108dbd7SJacob Faibussowitsch   PetscFunctionReturn(0);
25f108dbd7SJacob Faibussowitsch }
26c4762a1bSJed Brown 
27f108dbd7SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
28f108dbd7SJacob Faibussowitsch {
29f108dbd7SJacob Faibussowitsch   PetscFunctionBeginUser;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc, char **argv)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   DM             dm;
40f108dbd7SJacob Faibussowitsch   DMLabel        OQLabel;
41f108dbd7SJacob Faibussowitsch   Vec            OQ;
42c4762a1bSJed Brown   AppCtx         ctx;
43c4762a1bSJed Brown 
44*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&OQ));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
51*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
52*b122ec5aSJacob 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}}
61e600fa54SMatthew 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