xref: /petsc/src/dm/impls/plex/tests/ex32.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
1c4762a1bSJed Brown static char help[] = "Tests for periodic mesh output\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
CheckMesh(DM dm)5d71ae5a4SJacob Faibussowitsch PetscErrorCode CheckMesh(DM dm)
6d71ae5a4SJacob Faibussowitsch {
75f80ce2aSJacob Faibussowitsch   PetscReal detJ, J[9];
8c4762a1bSJed Brown   PetscReal vol;
95f80ce2aSJacob Faibussowitsch   PetscInt  dim, depth, cStart, cEnd, c;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
139566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
15c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
169566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, NULL, J, NULL, &detJ));
1763a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %" PetscInt_FMT " is inverted, |J| = %g", c, (double)detJ);
18c4762a1bSJed Brown     if (depth > 1) {
199566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL));
2063a3b9bcSJacob Faibussowitsch       PetscCheck(vol > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %" PetscInt_FMT " is inverted, vol = %g", c, (double)vol);
21c4762a1bSJed Brown     }
22c4762a1bSJed Brown   }
23*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,DM * dm)26d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
27d71ae5a4SJacob Faibussowitsch {
28c4762a1bSJed Brown   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
309566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
319566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
329566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
33*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
main(int argc,char ** argv)36d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
37d71ae5a4SJacob Faibussowitsch {
3830602db0SMatthew G. Knepley   DM dm;
39c4762a1bSJed Brown 
40327415f7SBarry Smith   PetscFunctionBeginUser;
419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
429566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
439566063dSJacob Faibussowitsch   PetscCall(CheckMesh(dm));
449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
459566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
46b122ec5aSJacob Faibussowitsch   return 0;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*TEST
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   test:
52c4762a1bSJed Brown     suffix: 0
5330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -dm_view ::ascii_info_detail
54c4762a1bSJed Brown   test:
55c4762a1bSJed Brown     suffix: 1
56c4762a1bSJed Brown     nsize: 2
57e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
58c4762a1bSJed Brown   test:
59c4762a1bSJed Brown     suffix: 2
60c4762a1bSJed Brown     nsize: 2
61e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
62c4762a1bSJed Brown   test:
63c4762a1bSJed Brown     suffix: 3
64c4762a1bSJed Brown     nsize: 4
65e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
66c4762a1bSJed Brown   test:
67c4762a1bSJed Brown     suffix: 4
68c4762a1bSJed Brown     nsize: 2
69e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,1,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
70c4762a1bSJed Brown   test:
71c4762a1bSJed Brown     suffix: 5
72c4762a1bSJed Brown     nsize: 2
73e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
741fe9daf4SMatthew G. Knepley   # This checks that the SF with extra root for periodic cut still checks
751fe9daf4SMatthew G. Knepley   test:
761fe9daf4SMatthew G. Knepley     suffix: 5_hdf5
771fe9daf4SMatthew G. Knepley     requires: hdf5
781fe9daf4SMatthew G. Knepley     nsize: 2
791fe9daf4SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view hdf5:mesh.h5
80c4762a1bSJed Brown   test:
81c4762a1bSJed Brown     suffix: 6
82c4762a1bSJed Brown     nsize: 4
83e600fa54SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 6,2,0 -dm_plex_box_bd periodic,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
84c4762a1bSJed Brown 
85c4762a1bSJed Brown TEST*/
86