1 static char help[] = "Element closure restrictions in tensor/lexicographic/spectral-element ordering using DMPlex\n\n"; 2 3 #include <petscdmplex.h> 4 5 int main(int argc, char **argv) 6 { 7 DM dm; 8 PetscSection section; 9 PetscFE fe; 10 PetscInt dim,c,cStart,cEnd; 11 PetscErrorCode ierr; 12 13 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 14 ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr); 15 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 16 ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 17 ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 18 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 19 20 ierr = PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe);CHKERRQ(ierr); 21 ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 22 ierr = DMAddField(dm,NULL,(PetscObject)fe);CHKERRQ(ierr); 23 ierr = DMCreateDS(dm);CHKERRQ(ierr); 24 ierr = DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL);CHKERRQ(ierr); 25 ierr = DMGetLocalSection(dm,§ion);CHKERRQ(ierr); 26 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 27 for (c=cStart; c<cEnd; c++) { 28 PetscInt numindices,*indices; 29 ierr = DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);CHKERRQ(ierr); 30 ierr = PetscPrintf(PETSC_COMM_SELF,"Element #%D\n",c-cStart);CHKERRQ(ierr); 31 ierr = PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 32 ierr = DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);CHKERRQ(ierr); 33 } 34 35 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 36 ierr = DMDestroy(&dm);CHKERRQ(ierr); 37 ierr = PetscFinalize(); 38 return ierr; 39 } 40 41 /*TEST 42 43 test: 44 suffix: 1d_q2 45 args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 46 test: 47 suffix: 2d_q1 48 args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 49 test: 50 suffix: 2d_q2 51 args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 52 test: 53 suffix: 2d_q3 54 args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 55 test: 56 suffix: 3d_q1 57 args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 58 59 TEST*/ 60