1*e7c9d70eSJeremy Theler static const char help[] = "Test for non-manifold interpolation"; 2*e7c9d70eSJeremy Theler 3*e7c9d70eSJeremy Theler #include <petscdmplex.h> 4*e7c9d70eSJeremy Theler 5*e7c9d70eSJeremy Theler /* 6*e7c9d70eSJeremy Theler 3-------------7 7*e7c9d70eSJeremy Theler /| /| 8*e7c9d70eSJeremy Theler / | / | 9*e7c9d70eSJeremy Theler / | / | 10*e7c9d70eSJeremy Theler 1-------------5 | 11*e7c9d70eSJeremy Theler | | | | 12*e7c9d70eSJeremy Theler | | | | 13*e7c9d70eSJeremy Theler | | | | 14*e7c9d70eSJeremy Theler | | | | 15*e7c9d70eSJeremy Theler z 4---------|---8 16*e7c9d70eSJeremy Theler ^ / | / 17*e7c9d70eSJeremy Theler | y | / 18*e7c9d70eSJeremy Theler |/ |/ 19*e7c9d70eSJeremy Theler 2--->-x-------6-------------9 20*e7c9d70eSJeremy Theler */ 21*e7c9d70eSJeremy Theler int main(int argc, char **argv) 22*e7c9d70eSJeremy Theler { 23*e7c9d70eSJeremy Theler DM dm, idm; 24*e7c9d70eSJeremy Theler DMLabel ctLabel; 25*e7c9d70eSJeremy Theler PetscBool has_vtk = PETSC_FALSE; 26*e7c9d70eSJeremy Theler 27*e7c9d70eSJeremy Theler // 9 vertices 28*e7c9d70eSJeremy Theler // 1 edge 29*e7c9d70eSJeremy Theler // 0 faces 30*e7c9d70eSJeremy Theler // 1 volume 31*e7c9d70eSJeremy Theler PetscInt num_points[4] = {9, 1, 0, 1}; 32*e7c9d70eSJeremy Theler 33*e7c9d70eSJeremy Theler // point 0 = hexahedron (defined by 8 vertices) 34*e7c9d70eSJeremy Theler // points 1-9 = vertices 35*e7c9d70eSJeremy Theler // point 10 = edged (defined by 2 vertices) 36*e7c9d70eSJeremy Theler PetscInt cone_size[11] = {8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2}; 37*e7c9d70eSJeremy Theler 38*e7c9d70eSJeremy Theler // hexahedron defined by points 39*e7c9d70eSJeremy Theler PetscInt cones[11] = {3, 4, 2, 1, 7, 5, 6, 8, 6, 9}; 40*e7c9d70eSJeremy Theler PetscInt cone_orientations[11] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 41*e7c9d70eSJeremy Theler PetscScalar vertex_coords[3 * 9] = {0, 0, 1, 0, 0, 0, 0, 1, 1, 0, 1, 0, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 1, 0, 2, 0, 0}; 42*e7c9d70eSJeremy Theler 43*e7c9d70eSJeremy Theler PetscFunctionBeginUser; 44*e7c9d70eSJeremy Theler PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 45*e7c9d70eSJeremy Theler PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Output VTK?", "ex66"); 46*e7c9d70eSJeremy Theler PetscCall(PetscOptionsGetBool(NULL, NULL, "-vtk", &has_vtk, NULL)); 47*e7c9d70eSJeremy Theler PetscOptionsEnd(); 48*e7c9d70eSJeremy Theler 49*e7c9d70eSJeremy Theler PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 50*e7c9d70eSJeremy Theler PetscCall(PetscObjectSetName((PetscObject)dm, "cubeline-fromdag")); 51*e7c9d70eSJeremy Theler PetscCall(DMSetType(dm, DMPLEX)); 52*e7c9d70eSJeremy Theler PetscCall(DMSetDimension(dm, 3)); 53*e7c9d70eSJeremy Theler 54*e7c9d70eSJeremy Theler PetscCall(DMPlexCreateFromDAG(dm, 3, num_points, cone_size, cones, cone_orientations, vertex_coords)); 55*e7c9d70eSJeremy Theler PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 56*e7c9d70eSJeremy Theler 57*e7c9d70eSJeremy Theler // TODO: make it work with a DM made from a msh file 58*e7c9d70eSJeremy Theler // PetscCall(DMPlexCreateGmshFromFile(PETSC_COMM_WORLD, "cube-line.msh", PETSC_FALSE, &dm)); 59*e7c9d70eSJeremy Theler // PetscCall(PetscObjectSetName((PetscObject)dm, "msh")); 60*e7c9d70eSJeremy Theler 61*e7c9d70eSJeremy Theler // Must set cell types 62*e7c9d70eSJeremy Theler PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel)); 63*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 0, DM_POLYTOPE_HEXAHEDRON)); 64*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 1, DM_POLYTOPE_POINT)); 65*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 2, DM_POLYTOPE_POINT)); 66*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 3, DM_POLYTOPE_POINT)); 67*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 4, DM_POLYTOPE_POINT)); 68*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 5, DM_POLYTOPE_POINT)); 69*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 6, DM_POLYTOPE_POINT)); 70*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 7, DM_POLYTOPE_POINT)); 71*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 8, DM_POLYTOPE_POINT)); 72*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 9, DM_POLYTOPE_POINT)); 73*e7c9d70eSJeremy Theler PetscCall(DMLabelSetValue(ctLabel, 10, DM_POLYTOPE_SEGMENT)); 74*e7c9d70eSJeremy Theler 75*e7c9d70eSJeremy Theler // interpolate (make sure to use -interp_dm_plex_stratify_celltype) 76*e7c9d70eSJeremy Theler PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "interp_")); 77*e7c9d70eSJeremy Theler PetscCall(DMPlexInterpolate(dm, &idm)); 78*e7c9d70eSJeremy Theler PetscCall(DMDestroy(&dm)); 79*e7c9d70eSJeremy Theler dm = idm; 80*e7c9d70eSJeremy Theler 81*e7c9d70eSJeremy Theler PetscCall(DMSetFromOptions(dm)); 82*e7c9d70eSJeremy Theler PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 83*e7c9d70eSJeremy Theler 84*e7c9d70eSJeremy Theler if (has_vtk) { 85*e7c9d70eSJeremy Theler PetscViewer viewer; 86*e7c9d70eSJeremy Theler PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 87*e7c9d70eSJeremy Theler PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK)); 88*e7c9d70eSJeremy Theler PetscCall(PetscViewerFileSetName(viewer, "ex66.vtk")); 89*e7c9d70eSJeremy Theler PetscCall(DMView(dm, viewer)); 90*e7c9d70eSJeremy Theler PetscCall(PetscViewerDestroy(&viewer)); 91*e7c9d70eSJeremy Theler } 92*e7c9d70eSJeremy Theler 93*e7c9d70eSJeremy Theler PetscCall(DMDestroy(&dm)); 94*e7c9d70eSJeremy Theler PetscCall(PetscFinalize()); 95*e7c9d70eSJeremy Theler return 0; 96*e7c9d70eSJeremy Theler } 97*e7c9d70eSJeremy Theler 98*e7c9d70eSJeremy Theler /*TEST 99*e7c9d70eSJeremy Theler test: 100*e7c9d70eSJeremy Theler suffix: 0 101*e7c9d70eSJeremy Theler args: -interp_dm_plex_stratify_celltype -dm_view ::ascii_info_detail -interp_dm_view ::ascii_info_detail 102*e7c9d70eSJeremy Theler 103*e7c9d70eSJeremy Theler TEST*/ 104