xref: /petsc/src/dm/impls/plex/tutorials/ex8.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
5*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewOffsets(DM dm, Vec X)
6*d71ae5a4SJacob Faibussowitsch {
7dd7309edSJed Brown   PetscInt           num_elem, elem_size, num_comp, num_dof;
8dd7309edSJed Brown   PetscInt          *elem_restr_offsets;
9dd7309edSJed Brown   const PetscScalar *x = NULL;
10dd7309edSJed Brown   const char        *name;
11dd7309edSJed Brown 
12dd7309edSJed Brown   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
159371c9d4SSatish 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));
169566063dSJacob Faibussowitsch   if (X) PetscCall(VecGetArrayRead(X, &x));
17dd7309edSJed Brown   for (PetscInt c = 0; c < num_elem; c++) {
189566063dSJacob Faibussowitsch     PetscCall(PetscIntView(elem_size, &elem_restr_offsets[c * elem_size], PETSC_VIEWER_STDOUT_SELF));
19dd7309edSJed Brown     if (x) {
2048a46eb9SPierre 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));
21dd7309edSJed Brown     }
22dd7309edSJed Brown   }
239566063dSJacob Faibussowitsch   if (X) PetscCall(VecRestoreArrayRead(X, &x));
249566063dSJacob Faibussowitsch   PetscCall(PetscFree(elem_restr_offsets));
25dd7309edSJed Brown   PetscFunctionReturn(0);
26dd7309edSJed Brown }
27dd7309edSJed Brown 
28*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
29*d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown   DM           dm;
31c4762a1bSJed Brown   PetscSection section;
32c4762a1bSJed Brown   PetscFE      fe;
3330602db0SMatthew G. Knepley   PetscInt     dim, c, cStart, cEnd;
348ba16be8SJed Brown   PetscBool    view_coord = PETSC_FALSE, tensor = PETSC_TRUE, project = PETSC_FALSE;
35c4762a1bSJed Brown 
36327415f7SBarry Smith   PetscFunctionBeginUser;
379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
38d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX");
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL));
408ba16be8SJed Brown   PetscCall(PetscOptionsBool("-project_coordinates", "Call DMProjectCoordinates explicitly", "ex8.c", project, &project, NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL));
42d0609cedSBarry Smith   PetscOptionsEnd();
4346181b2aSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
459566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
469566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
478ba16be8SJed Brown   if (project) {
488ba16be8SJed Brown     PetscFE  fe_coords;
498ba16be8SJed Brown     PetscInt cdim;
508ba16be8SJed Brown     PetscCall(DMGetCoordinateDim(dm, &cdim));
518ba16be8SJed Brown     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, cdim, cdim, PETSC_FALSE, 1, 1, &fe_coords));
528ba16be8SJed Brown     PetscCall(DMProjectCoordinates(dm, fe_coords));
538ba16be8SJed Brown     PetscCall(PetscFEDestroy(&fe_coords));
548ba16be8SJed Brown   }
559566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DETERMINE, &fe));
599566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
609566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
619566063dSJacob Faibussowitsch   if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
64c4762a1bSJed Brown   for (c = cStart; c < cEnd; c++) {
65c4762a1bSJed Brown     PetscInt numindices, *indices;
669566063dSJacob Faibussowitsch     PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL));
679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT "\n", c - cStart));
689566063dSJacob Faibussowitsch     PetscCall(PetscIntView(numindices, indices, PETSC_VIEWER_STDOUT_SELF));
699566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL));
70c4762a1bSJed Brown   }
7146181b2aSJed Brown   if (view_coord) {
7246181b2aSJed Brown     DM       cdm;
7346181b2aSJed Brown     Vec      X;
7446181b2aSJed Brown     PetscInt cdim;
758fb5bd83SMatthew G. Knepley 
768fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &cdim));
779566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
789566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)cdm, "coords"));
799566063dSJacob Faibussowitsch     if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
808fb5bd83SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
818fb5bd83SMatthew G. Knepley       const PetscScalar *array;
82dd7309edSJed Brown       PetscScalar       *x = NULL;
8346181b2aSJed Brown       PetscInt           ndof;
848fb5bd83SMatthew G. Knepley       PetscBool          isDG;
858fb5bd83SMatthew G. Knepley 
868fb5bd83SMatthew G. Knepley       PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x));
8746181b2aSJed Brown       PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim");
889566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c - cStart));
898fb5bd83SMatthew G. Knepley       for (PetscInt i = 0; i < ndof; i += cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF));
908fb5bd83SMatthew G. Knepley       PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x));
9146181b2aSJed Brown     }
929566063dSJacob Faibussowitsch     PetscCall(ViewOffsets(dm, NULL));
938fb5bd83SMatthew G. Knepley     PetscCall(DMGetCoordinatesLocal(dm, &X));
949566063dSJacob Faibussowitsch     PetscCall(ViewOffsets(cdm, X));
9546181b2aSJed Brown   }
969566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
979566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
99b122ec5aSJacob Faibussowitsch   return 0;
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /*TEST
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   test:
105c4762a1bSJed Brown     suffix: 1d_q2
10630602db0SMatthew G. Knepley     args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2
107c4762a1bSJed Brown   test:
108c4762a1bSJed Brown     suffix: 2d_q1
10930602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
110c4762a1bSJed Brown   test:
111c4762a1bSJed Brown     suffix: 2d_q2
11230602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2
113c4762a1bSJed Brown   test:
114c4762a1bSJed Brown     suffix: 2d_q3
11530602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1
116c4762a1bSJed Brown   test:
117c4762a1bSJed Brown     suffix: 3d_q1
11830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1
11946181b2aSJed Brown   test:
12046181b2aSJed Brown     suffix: 1d_q1_periodic
121a05c9aa3SJed Brown     requires: !complex
122dd7309edSJed 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
12346181b2aSJed Brown   test:
12446181b2aSJed Brown     suffix: 2d_q1_periodic
125a05c9aa3SJed Brown     requires: !complex
126dd7309edSJed 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
12746181b2aSJed Brown   test:
128a05c9aa3SJed Brown     suffix: 3d_q1_periodic
129a05c9aa3SJed Brown     requires: !complex
130a05c9aa3SJed 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
131a05c9aa3SJed Brown   test:
1328ba16be8SJed Brown     suffix: 3d_q1_periodic_project
1338ba16be8SJed Brown     requires: !complex
1348ba16be8SJed 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
1358ba16be8SJed Brown 
1368ba16be8SJed Brown   test:
137a05c9aa3SJed Brown     suffix: 3d_q2_periodic  # not actually periodic because only 2 cells
13846181b2aSJed 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
139c4762a1bSJed Brown 
140c4762a1bSJed Brown TEST*/
141