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