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 5dd7309edSJed Brown static PetscErrorCode ViewOffsets(DM dm, Vec X) 6dd7309edSJed Brown { 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)); 15d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF,"DM %s offsets: num_elem %" PetscInt_FMT ", size %" PetscInt_FMT 16d0609cedSBarry Smith ", 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) { 21dd7309edSJed Brown for (PetscInt i=0; i<elem_size; i++) { 229566063dSJacob Faibussowitsch PetscCall(PetscScalarView(num_comp, &x[elem_restr_offsets[c*elem_size+i]], PETSC_VIEWER_STDERR_SELF)); 23dd7309edSJed Brown } 24dd7309edSJed Brown } 25dd7309edSJed Brown } 269566063dSJacob Faibussowitsch if (X) PetscCall(VecRestoreArrayRead(X, &x)); 279566063dSJacob Faibussowitsch PetscCall(PetscFree(elem_restr_offsets)); 28dd7309edSJed Brown PetscFunctionReturn(0); 29dd7309edSJed Brown } 30dd7309edSJed Brown 31c4762a1bSJed Brown int main(int argc, char **argv) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown DM dm; 34c4762a1bSJed Brown PetscSection section; 35c4762a1bSJed Brown PetscFE fe; 3630602db0SMatthew G. Knepley PetscInt dim,c,cStart,cEnd; 37*8ba16be8SJed Brown PetscBool view_coord = PETSC_FALSE, tensor = PETSC_TRUE, project = PETSC_FALSE; 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 40d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Tensor closure restrictions", "DMPLEX"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-closure_tensor", "Apply DMPlexSetClosurePermutationTensor", "ex8.c", tensor, &tensor, NULL)); 42*8ba16be8SJed Brown PetscCall(PetscOptionsBool("-project_coordinates", "Call DMProjectCoordinates explicitly", "ex8.c", project, &project, NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL)); 44d0609cedSBarry Smith PetscOptionsEnd(); 4546181b2aSJed Brown 469566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 479566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 49*8ba16be8SJed Brown if (project) { 50*8ba16be8SJed Brown PetscFE fe_coords; 51*8ba16be8SJed Brown PetscInt cdim; 52*8ba16be8SJed Brown PetscCall(DMGetCoordinateDim(dm, &cdim)); 53*8ba16be8SJed Brown PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, cdim, cdim, PETSC_FALSE, 1, 1, &fe_coords)); 54*8ba16be8SJed Brown PetscCall(DMProjectCoordinates(dm, fe_coords)); 55*8ba16be8SJed Brown PetscCall(PetscFEDestroy(&fe_coords)); 56*8ba16be8SJed Brown } 579566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe)); 619566063dSJacob Faibussowitsch PetscCall(DMAddField(dm,NULL,(PetscObject)fe)); 629566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 639566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL)); 649566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm,§ion)); 659566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 66c4762a1bSJed Brown for (c=cStart; c<cEnd; c++) { 67c4762a1bSJed Brown PetscInt numindices,*indices; 689566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT "\n",c-cStart)); 709566063dSJacob Faibussowitsch PetscCall(PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF)); 719566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 72c4762a1bSJed Brown } 7346181b2aSJed Brown if (view_coord) { 7446181b2aSJed Brown DM cdm; 7546181b2aSJed Brown Vec X; 7646181b2aSJed Brown PetscInt cdim; 778fb5bd83SMatthew G. Knepley 788fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 799566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) cdm, "coords")); 819566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 828fb5bd83SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 838fb5bd83SMatthew G. Knepley const PetscScalar *array; 84dd7309edSJed Brown PetscScalar *x = NULL; 8546181b2aSJed Brown PetscInt ndof; 868fb5bd83SMatthew G. Knepley PetscBool isDG; 878fb5bd83SMatthew G. Knepley 888fb5bd83SMatthew G. Knepley PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 8946181b2aSJed Brown PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim"); 909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c-cStart)); 918fb5bd83SMatthew G. Knepley for (PetscInt i = 0; i < ndof; i+= cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF)); 928fb5bd83SMatthew G. Knepley PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 9346181b2aSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(ViewOffsets(dm, NULL)); 958fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &X)); 969566063dSJacob Faibussowitsch PetscCall(ViewOffsets(cdm, X)); 9746181b2aSJed Brown } 989566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 101b122ec5aSJacob Faibussowitsch return 0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /*TEST 105c4762a1bSJed Brown 106c4762a1bSJed Brown test: 107c4762a1bSJed Brown suffix: 1d_q2 10830602db0SMatthew G. Knepley args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 109c4762a1bSJed Brown test: 110c4762a1bSJed Brown suffix: 2d_q1 11130602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 112c4762a1bSJed Brown test: 113c4762a1bSJed Brown suffix: 2d_q2 11430602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 115c4762a1bSJed Brown test: 116c4762a1bSJed Brown suffix: 2d_q3 11730602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 118c4762a1bSJed Brown test: 119c4762a1bSJed Brown suffix: 3d_q1 12030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 12146181b2aSJed Brown test: 12246181b2aSJed Brown suffix: 1d_q1_periodic 123a05c9aa3SJed Brown requires: !complex 124dd7309edSJed 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 12546181b2aSJed Brown test: 12646181b2aSJed Brown suffix: 2d_q1_periodic 127a05c9aa3SJed Brown requires: !complex 128dd7309edSJed 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 12946181b2aSJed Brown test: 130a05c9aa3SJed Brown suffix: 3d_q1_periodic 131a05c9aa3SJed Brown requires: !complex 132a05c9aa3SJed 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 133a05c9aa3SJed Brown test: 134*8ba16be8SJed Brown suffix: 3d_q1_periodic_project 135*8ba16be8SJed Brown requires: !complex 136*8ba16be8SJed 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 137*8ba16be8SJed Brown 138*8ba16be8SJed Brown test: 139a05c9aa3SJed Brown suffix: 3d_q2_periodic # not actually periodic because only 2 cells 14046181b2aSJed 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 141c4762a1bSJed Brown 142c4762a1bSJed Brown TEST*/ 143