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