xref: /petsc/src/dm/impls/plex/tutorials/ex8.c (revision 78f1d13956708d00685f8d108946ce2901b5ba49)
1 static char help[] = "Element closure restrictions in tensor/lexicographic/spectral-element ordering using DMPlex\n\n";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode ViewOffsets(DM dm, Vec X)
6 {
7   PetscInt          num_elem, elem_size, num_comp, num_dof;
8   PetscInt          *elem_restr_offsets;
9   const PetscScalar *x = NULL;
10   const char        *name;
11 
12   PetscFunctionBegin;
13   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
14   PetscCall(DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
15   PetscCall(PetscPrintf(PETSC_COMM_SELF,"DM %s offsets: num_elem %" PetscInt_FMT ", size %" PetscInt_FMT
16                         ", comp %" PetscInt_FMT ", dof %" PetscInt_FMT "\n",name, num_elem, elem_size, num_comp, num_dof));
17   if (X) PetscCall(VecGetArrayRead(X, &x));
18   for (PetscInt c=0; c<num_elem; c++) {
19     PetscCall(PetscIntView(elem_size, &elem_restr_offsets[c*elem_size], PETSC_VIEWER_STDOUT_SELF));
20     if (x) {
21       for (PetscInt i=0; i<elem_size; i++) {
22         PetscCall(PetscScalarView(num_comp, &x[elem_restr_offsets[c*elem_size+i]], PETSC_VIEWER_STDERR_SELF));
23       }
24     }
25   }
26   if (X) PetscCall(VecRestoreArrayRead(X, &x));
27   PetscCall(PetscFree(elem_restr_offsets));
28   PetscFunctionReturn(0);
29 }
30 
31 int main(int argc, char **argv)
32 {
33   DM             dm;
34   PetscSection   section;
35   PetscFE        fe;
36   PetscInt       dim,c,cStart,cEnd;
37   PetscBool      view_coord = PETSC_FALSE, tensor = PETSC_TRUE;
38 
39   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
40   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX");
41   PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL));
42   PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL));
43   PetscOptionsEnd();
44 
45   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
46   PetscCall(DMSetType(dm, DMPLEX));
47   PetscCall(DMSetFromOptions(dm));
48   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
49   PetscCall(DMGetDimension(dm, &dim));
50 
51   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe));
52   PetscCall(DMAddField(dm,NULL,(PetscObject)fe));
53   PetscCall(DMCreateDS(dm));
54   if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL));
55   PetscCall(DMGetLocalSection(dm,&section));
56   PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
57   for (c=cStart; c<cEnd; c++) {
58     PetscInt numindices,*indices;
59     PetscCall(DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL));
60     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT "\n",c-cStart));
61     PetscCall(PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF));
62     PetscCall(DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL));
63   }
64   if (view_coord) {
65     DM cdm;
66     Vec X;
67     PetscInt cdim;
68     PetscCall(DMGetCoordinateDM(dm,&cdm));
69     PetscCall(PetscObjectSetName((PetscObject)cdm, "coords"));
70     if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm,PETSC_DETERMINE,NULL));
71     PetscCall(DMGetCoordinatesLocal(dm, &X));
72     PetscCall(DMGetCoordinateDim(dm, &cdim));
73     for (c=cStart; c<cEnd; c++) {
74       PetscScalar *x = NULL;
75       PetscInt ndof;
76       PetscCall(DMPlexVecGetClosure(cdm, NULL, X, c, &ndof, &x));
77       PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim");
78       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT " coordinates\n",c-cStart));
79       for (PetscInt i=0; i<ndof; i+= cdim) {
80         PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF));
81       }
82       PetscCall(DMPlexVecRestoreClosure(cdm, NULL, X, c, &ndof, &x));
83     }
84     PetscCall(ViewOffsets(dm, NULL));
85     PetscCall(ViewOffsets(cdm, X));
86   }
87   PetscCall(PetscFEDestroy(&fe));
88   PetscCall(DMDestroy(&dm));
89   PetscCall(PetscFinalize());
90   return 0;
91 }
92 
93 /*TEST
94 
95   test:
96     suffix: 1d_q2
97     args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2
98   test:
99     suffix: 2d_q1
100     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
101   test:
102     suffix: 2d_q2
103     args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
104   test:
105     suffix: 2d_q3
106     args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1
107   test:
108     suffix: 3d_q1
109     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
110   test:
111     suffix: 1d_q1_periodic
112     requires: !complex
113     args: -dm_plex_dim 1 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3 -dm_plex_box_bd periodic -dm_view -view_coord
114   test:
115     suffix: 2d_q1_periodic
116     requires: !complex
117     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3,2 -dm_plex_box_bd periodic,none -dm_view -view_coord
118   test:
119     suffix: 3d_q1_periodic
120     requires: !complex
121     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 3,2,1 -dm_plex_box_bd periodic,none,none -dm_view -view_coord
122   test:
123     suffix: 3d_q2_periodic  # not actually periodic because only 2 cells
124     args: -dm_plex_dim 3 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2 -dm_plex_box_bd periodic,none,periodic -dm_view
125 
126 TEST*/
127