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> 43c82e914SMatthew 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, §ion)); 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 } 843c82e914SMatthew G. Knepley { 853c82e914SMatthew G. Knepley DMLabel domain_label; 863c82e914SMatthew G. Knepley IS valueIS; 873c82e914SMatthew G. Knepley const PetscInt *values; 883c82e914SMatthew G. Knepley PetscInt Nv; 893c82e914SMatthew G. Knepley 903c82e914SMatthew G. Knepley PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 913c82e914SMatthew G. Knepley PetscCall(DMLabelGetNumValues(domain_label, &Nv)); 923c82e914SMatthew G. Knepley // Check that discretization bas any values on faces 933c82e914SMatthew G. Knepley { 943c82e914SMatthew G. Knepley PetscDS ds; 953c82e914SMatthew G. Knepley PetscFE fe; 963c82e914SMatthew G. Knepley IS fieldIS; 973c82e914SMatthew G. Knepley const PetscInt *fields; 983c82e914SMatthew G. Knepley PetscInt Nf, dmf = 0, dsf = -1; 993c82e914SMatthew G. Knepley 1003c82e914SMatthew G. Knepley PetscCall(DMGetRegionDS(dm, domain_label, &fieldIS, &ds, NULL)); 1013c82e914SMatthew G. Knepley // Translate DM field number to DS field number 1023c82e914SMatthew G. Knepley PetscCall(ISGetIndices(fieldIS, &fields)); 1033c82e914SMatthew G. Knepley PetscCall(ISGetSize(fieldIS, &Nf)); 1043c82e914SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) { 1053c82e914SMatthew G. Knepley if (dmf == fields[f]) { 1063c82e914SMatthew G. Knepley dsf = f; 1073c82e914SMatthew G. Knepley break; 1083c82e914SMatthew G. Knepley } 1093c82e914SMatthew G. Knepley } 1103c82e914SMatthew G. Knepley PetscCall(ISRestoreIndices(fieldIS, &fields)); 1113c82e914SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, dsf, (PetscObject *)&fe)); 1123c82e914SMatthew G. Knepley PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe)); 1133c82e914SMatthew G. Knepley if (!fe) Nv = 0; 1143c82e914SMatthew G. Knepley } 1153c82e914SMatthew G. Knepley PetscCall(DMLabelGetValueIS(domain_label, &valueIS)); 1163c82e914SMatthew G. Knepley PetscCall(ISGetIndices(valueIS, &values)); 1173c82e914SMatthew G. Knepley for (PetscInt v = 0; v < Nv; ++v) { 1183c82e914SMatthew G. Knepley PetscInt *elem_restr_offsets_face; 1193c82e914SMatthew G. Knepley PetscInt num_elem, elem_size, num_comp, num_dof; 1203c82e914SMatthew G. Knepley 1213c82e914SMatthew G. Knepley PetscCall(DMPlexGetLocalOffsets(dm, domain_label, values[v], 1, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_face)); 1223c82e914SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "========= Face restriction; marker: %" PetscInt_FMT " ========\n", values[v])); 1233c82e914SMatthew 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)); 1243c82e914SMatthew 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)); 1253c82e914SMatthew G. Knepley PetscCall(PetscFree(elem_restr_offsets_face)); 1263c82e914SMatthew G. Knepley } 1273c82e914SMatthew G. Knepley PetscCall(ISRestoreIndices(valueIS, &values)); 1283c82e914SMatthew G. Knepley PetscCall(ISDestroy(&valueIS)); 1293c82e914SMatthew G. Knepley } 13046181b2aSJed Brown if (view_coord) { 131*04385073SMatthew G. Knepley DM cdm, cell_dm; 13246181b2aSJed Brown Vec X; 13346181b2aSJed Brown PetscInt cdim; 134*04385073SMatthew G. Knepley PetscBool sparseLocalize; 1358fb5bd83SMatthew G. Knepley 1369f4ada15SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 1378fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)cdm, "coords")); 1409566063dSJacob Faibussowitsch if (tensor) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 141*04385073SMatthew G. Knepley PetscCall(DMLocalizeCoordinates(dm)); 142*04385073SMatthew G. Knepley PetscCall(DMGetSparseLocalize(dm, &sparseLocalize)); 1438fb5bd83SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1448fb5bd83SMatthew G. Knepley const PetscScalar *array; 145dd7309edSJed Brown PetscScalar *x = NULL; 14646181b2aSJed Brown PetscInt ndof; 1478fb5bd83SMatthew G. Knepley PetscBool isDG; 1488fb5bd83SMatthew G. Knepley 1498fb5bd83SMatthew G. Knepley PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 15046181b2aSJed Brown PetscCheck(ndof % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "ndof not divisible by cdim"); 1519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element #%" PetscInt_FMT " coordinates\n", c - cStart)); 1528fb5bd83SMatthew G. Knepley for (PetscInt i = 0; i < ndof; i += cdim) PetscCall(PetscScalarView(cdim, &x[i], PETSC_VIEWER_STDOUT_SELF)); 1538fb5bd83SMatthew G. Knepley PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &ndof, &array, &x)); 15446181b2aSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(ViewOffsets(dm, NULL)); 1568fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &X)); 1579566063dSJacob Faibussowitsch PetscCall(ViewOffsets(cdm, X)); 158*04385073SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cell_dm)); 159*04385073SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)cell_dm, "cell coords")); 160*04385073SMatthew G. Knepley if (cell_dm && !sparseLocalize) { 161*04385073SMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &X)); 162*04385073SMatthew G. Knepley PetscCall(ViewOffsets(cell_dm, X)); 163*04385073SMatthew G. Knepley } 16446181b2aSJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 167b122ec5aSJacob Faibussowitsch return 0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown /*TEST 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown suffix: 1d_q2 17430602db0SMatthew G. Knepley args: -dm_plex_dim 1 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2 175c4762a1bSJed Brown test: 176c4762a1bSJed Brown suffix: 2d_q1 17730602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 178c4762a1bSJed Brown test: 1793c82e914SMatthew G. Knepley suffix: 2d_q1d 180e327e467SRezgar Shakeri args: -dm_plex_dim 2 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -petscdualspace_lagrange_continuity 0 181e327e467SRezgar Shakeri test: 1823c82e914SMatthew G. Knepley suffix: 2d_q1_q1d 183e327e467SRezgar 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 184e327e467SRezgar Shakeri test: 185c4762a1bSJed Brown suffix: 2d_q2 18630602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 187c4762a1bSJed Brown test: 1883c82e914SMatthew G. Knepley suffix: 2d_q2_q1 1893c82e914SMatthew 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 1903c82e914SMatthew G. Knepley test: 1913c82e914SMatthew G. Knepley suffix: 2d_q2_p0 1923c82e914SMatthew 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 1933c82e914SMatthew G. Knepley test: 1943c82e914SMatthew G. Knepley suffix: 2d_q2_q1d 1953c82e914SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -num_fields 2 \ 1963c82e914SMatthew G. Knepley -f0_petscspace_degree 2 \ 1973c82e914SMatthew G. Knepley -f1_petscspace_degree 1 -f1_petscdualspace_lagrange_continuity 0 1983c82e914SMatthew G. Knepley test: 1993c82e914SMatthew G. Knepley suffix: 2d_q2_p1d 2003c82e914SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,1 -num_fields 2 \ 2013c82e914SMatthew G. Knepley -f0_petscspace_degree 2 \ 2023c82e914SMatthew G. Knepley -f1_petscspace_degree 1 -f1_petscspace_poly_tensor 0 -f1_petscdualspace_lagrange_continuity 0 2033c82e914SMatthew G. Knepley test: 204c4762a1bSJed Brown suffix: 2d_q3 20530602db0SMatthew G. Knepley args: -dm_plex_dim 2 -petscspace_degree 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 2063c82e914SMatthew G. Knepley 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: 3d_q1 20930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 21046181b2aSJed Brown test: 21146181b2aSJed Brown suffix: 1d_q1_periodic 212a05c9aa3SJed Brown requires: !complex 213dd7309edSJed 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 214*04385073SMatthew G. Knepley testset: 21546181b2aSJed Brown suffix: 2d_q1_periodic 216a05c9aa3SJed Brown requires: !complex 217*04385073SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 3,2 -dm_plex_box_bd periodic,none \ 218*04385073SMatthew G. Knepley -petscspace_degree 1 -dm_view -view_coord 219*04385073SMatthew G. Knepley 220*04385073SMatthew G. Knepley test: 221*04385073SMatthew G. Knepley test: 222*04385073SMatthew G. Knepley suffix: sparse 223*04385073SMatthew G. Knepley args: -dm_sparse_localize false -dm_localize 0 224*04385073SMatthew G. Knepley 22546181b2aSJed Brown test: 226a05c9aa3SJed Brown suffix: 3d_q1_periodic 227a05c9aa3SJed Brown requires: !complex 228a05c9aa3SJed 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 229a05c9aa3SJed Brown test: 2308ba16be8SJed Brown suffix: 3d_q1_periodic_project 2318ba16be8SJed Brown requires: !complex 2328ba16be8SJed 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 2338ba16be8SJed Brown 2348ba16be8SJed Brown test: 235a05c9aa3SJed Brown suffix: 3d_q2_periodic # not actually periodic because only 2 cells 23646181b2aSJed 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 237c4762a1bSJed Brown 238c4762a1bSJed Brown TEST*/ 239