Actual source code: ex8.c
petsc-3.14.6 2021-03-30
1: static char help[] = "Element closure restrictions in tensor/lexicographic/spectral-element ordering using DMPlex\n\n";
3: #include <petscdmplex.h>
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;
13: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
14: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Spectral/tensor element restrictions",NULL);
15: PetscOptionsRangeInt("-dim","Topological dimension",NULL,dim,&dim,NULL,1,3);
16: tmp = dim;
17: PetscOptionsIntArray("-cells","Number of cells per dimension",NULL,cells,&tmp,NULL);
18: PetscOptionsEnd();
20: DMPlexCreateBoxMesh(PETSC_COMM_WORLD,dim,PETSC_FALSE,cells,NULL,NULL,NULL,PETSC_TRUE,&dm);
21: PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe);
22: DMSetFromOptions(dm);
23: DMAddField(dm,NULL,(PetscObject)fe);
24: DMCreateDS(dm);
25: DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL);
26: DMGetLocalSection(dm,§ion);
27: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
28: for (c=cStart; c<cEnd; c++) {
29: PetscInt numindices,*indices;
30: DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);
31: PetscPrintf(PETSC_COMM_SELF,"Element #%D\n",c-cStart);
32: PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF);
33: DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL);
34: }
36: PetscFEDestroy(&fe);
37: DMDestroy(&dm);
38: PetscFinalize();
39: return ierr;
40: }
42: /*TEST
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
60: TEST*/