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, project = PETSC_FALSE; 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("-project_coordinates", "Call DMProjectCoordinates explicitly", "ex8.c", project, &project, NULL)); 43 PetscCall(PetscOptionsBool("-view_coord", "View coordinates of element closures", "ex8.c", view_coord, &view_coord, NULL)); 44 PetscOptionsEnd(); 45 46 PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 47 PetscCall(DMSetType(dm, DMPLEX)); 48 PetscCall(DMSetFromOptions(dm)); 49 if (project) { 50 PetscFE fe_coords; 51 PetscInt cdim; 52 PetscCall(DMGetCoordinateDim(dm, &cdim)); 53 PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, cdim, cdim, PETSC_FALSE, 1, 1, &fe_coords)); 54 PetscCall(DMProjectCoordinates(dm, fe_coords)); 55 PetscCall(PetscFEDestroy(&fe_coords)); 56 } 57 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 58 PetscCall(DMGetDimension(dm, &dim)); 59 60 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF,dim,1,PETSC_FALSE,NULL,PETSC_DETERMINE,&fe)); 61 PetscCall(DMAddField(dm,NULL,(PetscObject)fe)); 62 PetscCall(DMCreateDS(dm)); 63 if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL)); 64 PetscCall(DMGetLocalSection(dm,§ion)); 65 PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 66 for (c=cStart; c<cEnd; c++) { 67 PetscInt numindices,*indices; 68 PetscCall(DMPlexGetClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 69 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Element #%" PetscInt_FMT "\n",c-cStart)); 70 PetscCall(PetscIntView(numindices,indices,PETSC_VIEWER_STDOUT_SELF)); 71 PetscCall(DMPlexRestoreClosureIndices(dm,section,section,c,PETSC_TRUE,&numindices,&indices,NULL,NULL)); 72 } 73 if (view_coord) { 74 DM cdm; 75 Vec X; 76 PetscInt cdim; 77 78 PetscCall(DMGetCoordinateDim(dm, &cdim)); 79 PetscCall(DMGetCoordinateDM(dm, &cdm)); 80 PetscCall(PetscObjectSetName((PetscObject) cdm, "coords")); 81 if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 82 for (c = cStart; c < cEnd; ++c) { 83 const PetscScalar *array; 84 PetscScalar *x = NULL; 85 PetscInt ndof; 86 PetscBool isDG; 87 88 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 89 PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim"); 90 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c-cStart)); 91 for (PetscInt i = 0; i < ndof; i+= cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF)); 92 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 93 } 94 PetscCall(ViewOffsets(dm, NULL)); 95 PetscCall(DMGetCoordinatesLocal(dm, &X)); 96 PetscCall(ViewOffsets(cdm, X)); 97 } 98 PetscCall(PetscFEDestroy(&fe)); 99 PetscCall(DMDestroy(&dm)); 100 PetscCall(PetscFinalize()); 101 return 0; 102 } 103 104 /*TEST 105 106 test: 107 suffix: 1d_q2 108 args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 109 test: 110 suffix: 2d_q1 111 args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 112 test: 113 suffix: 2d_q2 114 args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 115 test: 116 suffix: 2d_q3 117 args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 118 test: 119 suffix: 3d_q1 120 args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 121 test: 122 suffix: 1d_q1_periodic 123 requires: !complex 124 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 125 test: 126 suffix: 2d_q1_periodic 127 requires: !complex 128 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 129 test: 130 suffix: 3d_q1_periodic 131 requires: !complex 132 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 133 test: 134 suffix: 3d_q1_periodic_project 135 requires: !complex 136 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 138 test: 139 suffix: 3d_q2_periodic # not actually periodic because only 2 cells 140 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 141 142 TEST*/ 143