xref: /petsc/src/dm/impls/plex/tutorials/ex8.c (revision 3c82e914a54a3bac717edd0eae8b48049be382c2)
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>
4*3c82e914SMatthew G. Knepley #include <petscds.h>
5c4762a1bSJed Brown 
6d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewOffsets(DM dm, Vec X)
7d71ae5a4SJacob Faibussowitsch {
8dd7309edSJed Brown   PetscInt           num_elem, elem_size, num_comp, num_dof;
9dd7309edSJed Brown   PetscInt          *elem_restr_offsets;
10dd7309edSJed Brown   const PetscScalar *x = NULL;
11dd7309edSJed Brown   const char        *name;
12dd7309edSJed Brown 
13dd7309edSJed Brown   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
169371c9d4SSatish Balay   PetscCall(PetscPrintf(PETSC_COMM_SELF, "DM %s offsets: num_elem %" PetscInt_FMT ", size %" PetscInt_FMT ", comp %" PetscInt_FMT ", dof %" PetscInt_FMT "\n", name, num_elem, elem_size, num_comp, num_dof));
179566063dSJacob Faibussowitsch   if (X) PetscCall(VecGetArrayRead(X, &x));
18dd7309edSJed Brown   for (PetscInt c = 0; c < num_elem; c++) {
199566063dSJacob Faibussowitsch     PetscCall(PetscIntView(elem_size, &elem_restr_offsets[c * elem_size], PETSC_VIEWER_STDOUT_SELF));
20dd7309edSJed Brown     if (x) {
2148a46eb9SPierre Jolivet       for (PetscInt i = 0; i < elem_size; i++) PetscCall(PetscScalarView(num_comp, &x[elem_restr_offsets[c * elem_size + i]], PETSC_VIEWER_STDERR_SELF));
22dd7309edSJed Brown     }
23dd7309edSJed Brown   }
249566063dSJacob Faibussowitsch   if (X) PetscCall(VecRestoreArrayRead(X, &x));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(elem_restr_offsets));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27dd7309edSJed Brown }
28dd7309edSJed Brown 
29d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   DM           dm;
32c4762a1bSJed Brown   PetscSection section;
33c4762a1bSJed Brown   PetscFE      fe;
34e327e467SRezgar Shakeri   PetscInt     dim, Nf = 1, c, cStart, cEnd;
358ba16be8SJed Brown   PetscBool    view_coord = PETSC_FALSE, tensor = PETSC_TRUE, project = PETSC_FALSE;
36c4762a1bSJed Brown 
37327415f7SBarry Smith   PetscFunctionBeginUser;
389566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
39d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX");
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL));
41e44f6aebSMatthew G. Knepley   PetscCall(PetscOptionsBool("-project_coordinates", "Call DMSetCoordinateDisc() explicitly", "ex8.c", project, &project, NULL));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL));
43e327e467SRezgar Shakeri   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of fields to use", "ex8.c", Nf, &Nf, NULL, 1));
44d0609cedSBarry Smith   PetscOptionsEnd();
4546181b2aSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
479566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
489566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
498ba16be8SJed Brown   if (project) {
508ba16be8SJed Brown     PetscFE  fe_coords;
518ba16be8SJed Brown     PetscInt cdim;
528ba16be8SJed Brown     PetscCall(DMGetCoordinateDim(dm, &cdim));
538ba16be8SJed Brown     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, cdim, cdim, PETSC_FALSE, 1, 1, &fe_coords));
54e44f6aebSMatthew G. Knepley     PetscCall(DMSetCoordinateDisc(dm, fe_coords, PETSC_TRUE));
558ba16be8SJed Brown     PetscCall(PetscFEDestroy(&fe_coords));
568ba16be8SJed Brown   }
579566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
589566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
59c4762a1bSJed Brown 
60e327e467SRezgar Shakeri   if (Nf == 1) {
619566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DETERMINE, &fe));
629566063dSJacob Faibussowitsch     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
63e327e467SRezgar Shakeri     PetscCall(PetscFEDestroy(&fe));
64e327e467SRezgar Shakeri   } else {
65e327e467SRezgar Shakeri     for (PetscInt f = 0; f < Nf; ++f) {
66e327e467SRezgar Shakeri       char prefix[16];
67e327e467SRezgar Shakeri       PetscCall(PetscSNPrintf(prefix, 16, "f%" PetscInt_FMT "_", f));
68e327e467SRezgar Shakeri       PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, prefix, PETSC_DETERMINE, &fe));
69e327e467SRezgar Shakeri       PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
70e327e467SRezgar Shakeri       PetscCall(PetscFEDestroy(&fe));
71e327e467SRezgar Shakeri     }
72e327e467SRezgar Shakeri   }
739566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
749566063dSJacob Faibussowitsch   if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
759566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
77c4762a1bSJed Brown   for (c = cStart; c < cEnd; c++) {
78c4762a1bSJed Brown     PetscInt numindices, *indices;
799566063dSJacob Faibussowitsch     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL));
809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT "\n", c - cStart));
819566063dSJacob Faibussowitsch     PetscCall(PetscIntView(numindices, indices, PETSC_VIEWER_STDOUT_SELF));
829566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL));
83c4762a1bSJed Brown   }
84*3c82e914SMatthew G. Knepley   {
85*3c82e914SMatthew G. Knepley     DMLabel         domain_label;
86*3c82e914SMatthew G. Knepley     IS              valueIS;
87*3c82e914SMatthew G. Knepley     const PetscInt *values;
88*3c82e914SMatthew G. Knepley     PetscInt        Nv;
89*3c82e914SMatthew G. Knepley 
90*3c82e914SMatthew G. Knepley     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
91*3c82e914SMatthew G. Knepley     PetscCall(DMLabelGetNumValues(domain_label, &Nv));
92*3c82e914SMatthew G. Knepley     // Check that discretization bas any values on faces
93*3c82e914SMatthew G. Knepley     {
94*3c82e914SMatthew G. Knepley       PetscDS         ds;
95*3c82e914SMatthew G. Knepley       PetscFE         fe;
96*3c82e914SMatthew G. Knepley       IS              fieldIS;
97*3c82e914SMatthew G. Knepley       const PetscInt *fields;
98*3c82e914SMatthew G. Knepley       PetscInt        Nf, dmf = 0, dsf = -1;
99*3c82e914SMatthew G. Knepley 
100*3c82e914SMatthew G. Knepley       PetscCall(DMGetRegionDS(dm, domain_label, &fieldIS, &ds, NULL));
101*3c82e914SMatthew G. Knepley       // Translate DM field number to DS field number
102*3c82e914SMatthew G. Knepley       PetscCall(ISGetIndices(fieldIS, &fields));
103*3c82e914SMatthew G. Knepley       PetscCall(ISGetSize(fieldIS, &Nf));
104*3c82e914SMatthew G. Knepley       for (PetscInt f = 0; f < Nf; ++f) {
105*3c82e914SMatthew G. Knepley         if (dmf == fields[f]) {
106*3c82e914SMatthew G. Knepley           dsf = f;
107*3c82e914SMatthew G. Knepley           break;
108*3c82e914SMatthew G. Knepley         }
109*3c82e914SMatthew G. Knepley       }
110*3c82e914SMatthew G. Knepley       PetscCall(ISRestoreIndices(fieldIS, &fields));
111*3c82e914SMatthew G. Knepley       PetscCall(PetscDSGetDiscretization(ds, dsf, (PetscObject *)&fe));
112*3c82e914SMatthew G. Knepley       PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe));
113*3c82e914SMatthew G. Knepley       if (!fe) Nv = 0;
114*3c82e914SMatthew G. Knepley     }
115*3c82e914SMatthew G. Knepley     PetscCall(DMLabelGetValueIS(domain_label, &valueIS));
116*3c82e914SMatthew G. Knepley     PetscCall(ISGetIndices(valueIS, &values));
117*3c82e914SMatthew G. Knepley     for (PetscInt v = 0; v < Nv; ++v) {
118*3c82e914SMatthew G. Knepley       PetscInt *elem_restr_offsets_face;
119*3c82e914SMatthew G. Knepley       PetscInt  num_elem, elem_size, num_comp, num_dof;
120*3c82e914SMatthew G. Knepley 
121*3c82e914SMatthew G. Knepley       PetscCall(DMPlexGetLocalOffsets(dm, domain_label, values[v], 1, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_face));
122*3c82e914SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "========= Face restriction; marker: %" PetscInt_FMT " ========\n", values[v]));
123*3c82e914SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "num_elem: %" PetscInt_FMT ", elem_size: %" PetscInt_FMT ", num_dof:  %" PetscInt_FMT "\n", num_elem, elem_size, num_dof));
124*3c82e914SMatthew G. Knepley       for (PetscInt c = 0; c < num_elem; c++) PetscCall(PetscIntView(elem_size, &elem_restr_offsets_face[c * elem_size], PETSC_VIEWER_STDOUT_SELF));
125*3c82e914SMatthew G. Knepley       PetscCall(PetscFree(elem_restr_offsets_face));
126*3c82e914SMatthew G. Knepley     }
127*3c82e914SMatthew G. Knepley     PetscCall(ISRestoreIndices(valueIS, &values));
128*3c82e914SMatthew G. Knepley     PetscCall(ISDestroy(&valueIS));
129*3c82e914SMatthew G. Knepley   }
13046181b2aSJed Brown   if (view_coord) {
13146181b2aSJed Brown     DM       cdm;
13246181b2aSJed Brown     Vec      X;
13346181b2aSJed Brown     PetscInt cdim;
1348fb5bd83SMatthew G. Knepley 
1359f4ada15SMatthew G. Knepley     PetscCall(DMGetCoordinatesLocalSetUp(dm));
1368fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &cdim));
1379566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
1389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)cdm, "coords"));
1399566063dSJacob Faibussowitsch     if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
1408fb5bd83SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
1418fb5bd83SMatthew G. Knepley       const PetscScalar *array;
142dd7309edSJed Brown       PetscScalar       *x = NULL;
14346181b2aSJed Brown       PetscInt           ndof;
1448fb5bd83SMatthew G. Knepley       PetscBool          isDG;
1458fb5bd83SMatthew G. Knepley 
1468fb5bd83SMatthew G. Knepley       PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x));
14746181b2aSJed Brown       PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim");
1489566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c - cStart));
1498fb5bd83SMatthew G. Knepley       for (PetscInt i = 0; i < ndof; i += cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF));
1508fb5bd83SMatthew G. Knepley       PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x));
15146181b2aSJed Brown     }
1529566063dSJacob Faibussowitsch     PetscCall(ViewOffsets(dm, NULL));
1538fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinatesLocal(dm, &X));
1549566063dSJacob Faibussowitsch     PetscCall(ViewOffsets(cdm, X));
15546181b2aSJed Brown   }
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
158b122ec5aSJacob Faibussowitsch   return 0;
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /*TEST
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   test:
164c4762a1bSJed Brown     suffix: 1d_q2
16530602db0SMatthew G. Knepley     args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2
166c4762a1bSJed Brown   test:
167c4762a1bSJed Brown     suffix: 2d_q1
16830602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
169c4762a1bSJed Brown   test:
170*3c82e914SMatthew G. Knepley     suffix: 2d_q1d
171e327e467SRezgar Shakeri     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -petscdualspace_lagrange_continuity 0
172e327e467SRezgar Shakeri   test:
173*3c82e914SMatthew G. Knepley     suffix: 2d_q1_q1d
174e327e467SRezgar Shakeri     args: -dm_plex_dim 2 -num_fields 2 -f0_petscspace_degree 1 -f1_petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -f1_petscdualspace_lagrange_continuity 0
175e327e467SRezgar Shakeri   test:
176c4762a1bSJed Brown     suffix: 2d_q2
17730602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
178c4762a1bSJed Brown   test:
179*3c82e914SMatthew G. Knepley     suffix: 2d_q2_q1
180*3c82e914SMatthew G. Knepley     args: -dm_plex_dim 2 -num_fields 2 -f0_petscspace_degree 2 -f1_petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,1
181*3c82e914SMatthew G. Knepley   test:
182*3c82e914SMatthew G. Knepley     suffix: 2d_q2_p0
183*3c82e914SMatthew G. Knepley     args: -dm_plex_dim 2 -num_fields 2 -f0_petscspace_degree 2 -f1_petscspace_degree 0 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -f1_petscdualspace_lagrange_continuity 0
184*3c82e914SMatthew G. Knepley   test:
185*3c82e914SMatthew G. Knepley     suffix: 2d_q2_q1d
186*3c82e914SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -num_fields 2 \
187*3c82e914SMatthew G. Knepley             -f0_petscspace_degree 2 \
188*3c82e914SMatthew G. Knepley             -f1_petscspace_degree 1 -f1_petscdualspace_lagrange_continuity 0
189*3c82e914SMatthew G. Knepley   test:
190*3c82e914SMatthew G. Knepley     suffix: 2d_q2_p1d
191*3c82e914SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -num_fields 2 \
192*3c82e914SMatthew G. Knepley             -f0_petscspace_degree 2 \
193*3c82e914SMatthew G. Knepley             -f1_petscspace_degree 1 -f1_petscspace_poly_tensor 0 -f1_petscdualspace_lagrange_continuity 0
194*3c82e914SMatthew G. Knepley   test:
195c4762a1bSJed Brown     suffix: 2d_q3
19630602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1
197*3c82e914SMatthew G. Knepley 
198c4762a1bSJed Brown   test:
199c4762a1bSJed Brown     suffix: 3d_q1
20030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
20146181b2aSJed Brown   test:
20246181b2aSJed Brown     suffix: 1d_q1_periodic
203a05c9aa3SJed Brown     requires: !complex
204dd7309edSJed Brown     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
20546181b2aSJed Brown   test:
20646181b2aSJed Brown     suffix: 2d_q1_periodic
207a05c9aa3SJed Brown     requires: !complex
208dd7309edSJed Brown     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
20946181b2aSJed Brown   test:
210a05c9aa3SJed Brown     suffix: 3d_q1_periodic
211a05c9aa3SJed Brown     requires: !complex
212a05c9aa3SJed Brown     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
213a05c9aa3SJed Brown   test:
2148ba16be8SJed Brown     suffix: 3d_q1_periodic_project
2158ba16be8SJed Brown     requires: !complex
2168ba16be8SJed Brown     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,3 -dm_plex_box_bd none,none,periodic -dm_view -view_coord -project_coordinates
2178ba16be8SJed Brown 
2188ba16be8SJed Brown   test:
219a05c9aa3SJed Brown     suffix: 3d_q2_periodic # not actually periodic because only 2 cells
22046181b2aSJed Brown     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
221c4762a1bSJed Brown 
222c4762a1bSJed Brown TEST*/
223