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 59371c9d4SSatish Balay static PetscErrorCode ViewOffsets(DM dm, Vec X) { 6dd7309edSJed Brown PetscInt num_elem, elem_size, num_comp, num_dof; 7dd7309edSJed Brown PetscInt *elem_restr_offsets; 8dd7309edSJed Brown const PetscScalar *x = NULL; 9dd7309edSJed Brown const char *name; 10dd7309edSJed Brown 11dd7309edSJed Brown PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 139566063dSJacob Faibussowitsch PetscCall(DMPlexGetLocalOffsets(dm, NULL, 0, 0, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets)); 149371c9d4SSatish 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)); 159566063dSJacob Faibussowitsch if (X) PetscCall(VecGetArrayRead(X, &x)); 16dd7309edSJed Brown for (PetscInt c = 0; c < num_elem; c++) { 179566063dSJacob Faibussowitsch PetscCall(PetscIntView(elem_size, &elem_restr_offsets[c * elem_size], PETSC_VIEWER_STDOUT_SELF)); 18dd7309edSJed Brown if (x) { 19*48a46eb9SPierre 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)); 20dd7309edSJed Brown } 21dd7309edSJed Brown } 229566063dSJacob Faibussowitsch if (X) PetscCall(VecRestoreArrayRead(X, &x)); 239566063dSJacob Faibussowitsch PetscCall(PetscFree(elem_restr_offsets)); 24dd7309edSJed Brown PetscFunctionReturn(0); 25dd7309edSJed Brown } 26dd7309edSJed Brown 279371c9d4SSatish Balay int main(int argc, char **argv) { 28c4762a1bSJed Brown DM dm; 29c4762a1bSJed Brown PetscSection section; 30c4762a1bSJed Brown PetscFE fe; 3130602db0SMatthew G. Knepley PetscInt dim, c, cStart, cEnd; 328ba16be8SJed Brown PetscBool view_coord = PETSC_FALSE, tensor = PETSC_TRUE, project = PETSC_FALSE; 33c4762a1bSJed Brown 34327415f7SBarry Smith PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 36d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX"); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL)); 388ba16be8SJed Brown PetscCall(PetscOptionsBool("-project_coordinates", "Call DMProjectCoordinates explicitly", "ex8.c", project, &project, NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL)); 40d0609cedSBarry Smith PetscOptionsEnd(); 4146181b2aSJed Brown 429566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 439566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 449566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 458ba16be8SJed Brown if (project) { 468ba16be8SJed Brown PetscFE fe_coords; 478ba16be8SJed Brown PetscInt cdim; 488ba16be8SJed Brown PetscCall(DMGetCoordinateDim(dm, &cdim)); 498ba16be8SJed Brown PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, cdim, cdim, PETSC_FALSE, 1, 1, &fe_coords)); 508ba16be8SJed Brown PetscCall(DMProjectCoordinates(dm, fe_coords)); 518ba16be8SJed Brown PetscCall(PetscFEDestroy(&fe_coords)); 528ba16be8SJed Brown } 539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 549566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DETERMINE, &fe)); 579566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 589566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 599566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 609566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 619566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 62c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) { 63c4762a1bSJed Brown PetscInt numindices, *indices; 649566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT "\n", c - cStart)); 669566063dSJacob Faibussowitsch PetscCall(PetscIntView(numindices, indices, PETSC_VIEWER_STDOUT_SELF)); 679566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, &numindices, &indices, NULL, NULL)); 68c4762a1bSJed Brown } 6946181b2aSJed Brown if (view_coord) { 7046181b2aSJed Brown DM cdm; 7146181b2aSJed Brown Vec X; 7246181b2aSJed Brown PetscInt cdim; 738fb5bd83SMatthew G. Knepley 748fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)cdm, "coords")); 779566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 788fb5bd83SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 798fb5bd83SMatthew G. Knepley const PetscScalar *array; 80dd7309edSJed Brown PetscScalar *x = NULL; 8146181b2aSJed Brown PetscInt ndof; 828fb5bd83SMatthew G. Knepley PetscBool isDG; 838fb5bd83SMatthew G. Knepley 848fb5bd83SMatthew G. Knepley PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 8546181b2aSJed Brown PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim"); 869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c - cStart)); 878fb5bd83SMatthew G. Knepley for (PetscInt i = 0; i < ndof; i += cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF)); 888fb5bd83SMatthew G. Knepley PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 8946181b2aSJed Brown } 909566063dSJacob Faibussowitsch PetscCall(ViewOffsets(dm, NULL)); 918fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &X)); 929566063dSJacob Faibussowitsch PetscCall(ViewOffsets(cdm, X)); 9346181b2aSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 97b122ec5aSJacob Faibussowitsch return 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /*TEST 101c4762a1bSJed Brown 102c4762a1bSJed Brown test: 103c4762a1bSJed Brown suffix: 1d_q2 10430602db0SMatthew G. Knepley args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 105c4762a1bSJed Brown test: 106c4762a1bSJed Brown suffix: 2d_q1 10730602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 108c4762a1bSJed Brown test: 109c4762a1bSJed Brown suffix: 2d_q2 11030602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 111c4762a1bSJed Brown test: 112c4762a1bSJed Brown suffix: 2d_q3 11330602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 114c4762a1bSJed Brown test: 115c4762a1bSJed Brown suffix: 3d_q1 11630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 11746181b2aSJed Brown test: 11846181b2aSJed Brown suffix: 1d_q1_periodic 119a05c9aa3SJed Brown requires: !complex 120dd7309edSJed 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 12146181b2aSJed Brown test: 12246181b2aSJed Brown suffix: 2d_q1_periodic 123a05c9aa3SJed Brown requires: !complex 124dd7309edSJed 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 12546181b2aSJed Brown test: 126a05c9aa3SJed Brown suffix: 3d_q1_periodic 127a05c9aa3SJed Brown requires: !complex 128a05c9aa3SJed 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 129a05c9aa3SJed Brown test: 1308ba16be8SJed Brown suffix: 3d_q1_periodic_project 1318ba16be8SJed Brown requires: !complex 1328ba16be8SJed 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 1338ba16be8SJed Brown 1348ba16be8SJed Brown test: 135a05c9aa3SJed Brown suffix: 3d_q2_periodic # not actually periodic because only 2 cells 13646181b2aSJed 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 137c4762a1bSJed Brown 138c4762a1bSJed Brown TEST*/ 139