xref: /petsc/src/dm/impls/plex/tests/ex27.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 static char help[] = "Test section ordering for FEM discretizations\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscds.h>
5 
6 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) {
7   PetscFunctionBegin;
8   PetscCall(DMCreate(comm, dm));
9   PetscCall(DMSetType(*dm, DMPLEX));
10   PetscCall(DMSetFromOptions(*dm));
11   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
12   PetscFunctionReturn(0);
13 }
14 
15 static PetscErrorCode TestLocalDofOrder(DM dm) {
16   PetscFE      fe[3];
17   PetscSection s;
18   PetscBool    simplex;
19   PetscInt     dim, Nf, f;
20 
21   PetscFunctionBegin;
22   PetscCall(DMGetDimension(dm, &dim));
23   PetscCall(DMPlexIsSimplex(dm, &simplex));
24   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field0_", -1, &fe[0]));
25   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "field1_", -1, &fe[1]));
26   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "field2_", -1, &fe[2]));
27 
28   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
29   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
30   PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe[2]));
31   PetscCall(DMCreateDS(dm));
32   PetscCall(DMGetLocalSection(dm, &s));
33   PetscCall(PetscObjectViewFromOptions((PetscObject)s, NULL, "-dof_view"));
34 
35   PetscCall(DMGetNumFields(dm, &Nf));
36   for (f = 0; f < Nf; ++f) PetscCall(PetscFEDestroy(&fe[f]));
37   PetscFunctionReturn(0);
38 }
39 
40 int main(int argc, char **argv) {
41   DM dm;
42 
43   PetscFunctionBeginUser;
44   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
45   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
46   PetscCall(TestLocalDofOrder(dm));
47   PetscCall(DMDestroy(&dm));
48   PetscCall(PetscFinalize());
49   return 0;
50 }
51 
52 /*TEST
53 
54   test:
55     suffix: tri_pm
56     requires: triangle
57     args: -dm_plex_box_faces 1,1 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view
58 
59   test:
60     suffix: quad_pm
61     requires:
62     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view
63 
64   test:
65     suffix: tri_fm
66     requires: triangle
67     args: -dm_coord_space 0 -dm_plex_box_faces 1,1 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view
68 
69   test:
70     suffix: quad_fm
71     requires:
72     args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view
73 
74 TEST*/
75