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