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