1c4762a1bSJed Brown static char help[] = "Tests interpolation and output of hybrid meshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscsf.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /* Much can be learned using 7c4762a1bSJed Brown -rd_dm_view -rd2_dm_view -rd_section_view -rd_vec_view -rd2_section_view */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown static PetscErrorCode redistribute_vec(DM dist_dm, PetscSF sf, Vec *v) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown DM dm, dist_v_dm; 12c4762a1bSJed Brown PetscSection section, dist_section; 13c4762a1bSJed Brown Vec dist_v; 14c4762a1bSJed Brown PetscMPIInt rank, size, p; 15c4762a1bSJed Brown 16c4762a1bSJed Brown PetscFunctionBegin; 175f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(*v, &dm)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(dm, §ion)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-rd_dm_view")); 205f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dist_dm, NULL, "-rd2_dm_view")); 21c4762a1bSJed Brown 225f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &dist_v_dm)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PetscObjectComm((PetscObject) *v), &dist_v)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetDM(dist_v, dist_v_dm)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) *v), &dist_section)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(dist_v_dm, dist_section)); 27c4762a1bSJed Brown 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) section, NULL, "-rd_section_view")); 295f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank)); 305f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 31c4762a1bSJed Brown for (p = 0; p < size; ++p) { 32c4762a1bSJed Brown if (p == rank) { 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) *v, NULL, "-rd_vec_view"));} 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBarrier((PetscObject) dm)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 36c4762a1bSJed Brown } 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeField(dm, sf, section, *v, dist_section, dist_v)); 38c4762a1bSJed Brown for (p = 0; p < size; ++p) { 39c4762a1bSJed Brown if (p == rank) { 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) dist_section, NULL, "-rd2_section_view")); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) dist_v, NULL, "-rd2_vec_view")); 42c4762a1bSJed Brown } 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBarrier((PetscObject) dm)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&dist_section)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dist_v_dm)); 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(v)); 51c4762a1bSJed Brown *v = dist_v; 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown static PetscErrorCode dm_view_geometry(DM dm, Vec cell_geom, Vec face_geom) 56c4762a1bSJed Brown { 57c4762a1bSJed Brown DM cell_dm, face_dm; 58c4762a1bSJed Brown PetscSection cell_section, face_section; 59c4762a1bSJed Brown const PetscScalar *cell_array, *face_array; 60c4762a1bSJed Brown const PetscInt *cells; 61c4762a1bSJed Brown PetscInt c, start_cell, end_cell; 62c4762a1bSJed Brown PetscInt f, start_face, end_face; 63c4762a1bSJed Brown PetscInt supportSize, offset; 64c4762a1bSJed Brown PetscMPIInt rank; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscFunctionBegin; 675f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* cells */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &start_cell, &end_cell)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(cell_geom, &cell_dm)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(cell_dm, &cell_section)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(cell_geom, &cell_array)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown for (c = start_cell; c < end_cell; ++c) { 76c4762a1bSJed Brown const PetscFVCellGeom *geom; 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(cell_section, c, &offset)); 78c4762a1bSJed Brown geom = (PetscFVCellGeom*)&cell_array[offset]; 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "rank %d c %D centroid %g,%g,%g vol %g\n", rank, c, (double)geom->centroid[0], (double)geom->centroid[1], (double)geom->centroid[2], (double)geom->volume)); 80c4762a1bSJed Brown } 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(cell_geom, &cell_array)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* faces */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 1, &start_face, &end_face)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(face_geom, &face_dm)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(face_dm, &face_section)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(face_geom, &face_array)); 89c4762a1bSJed Brown for (f = start_face; f < end_face; ++f) { 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupport(dm, f, &cells)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSupportSize(dm, f, &supportSize)); 92c4762a1bSJed Brown if (supportSize > 1) { 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(face_section, f, &offset)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "rank %d f %D cells %D,%D normal %g,%g,%g centroid %g,%g,%g\n", rank, f, cells[0], cells[1], (double) PetscRealPart(face_array[offset+0]), (double) PetscRealPart(face_array[offset+1]), (double) PetscRealPart(face_array[offset+2]), (double) PetscRealPart(face_array[offset+3]), (double) PetscRealPart(face_array[offset+4]), (double) PetscRealPart(face_array[offset+5]))); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown } 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(cell_geom, &cell_array)); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown int main(int argc, char **argv) 103c4762a1bSJed Brown { 10430602db0SMatthew G. Knepley DM dm, redist_dm; 105c4762a1bSJed Brown PetscPartitioner part; 10630602db0SMatthew G. Knepley PetscSF redist_sf; 107c4762a1bSJed Brown Vec cell_geom, face_geom; 10830602db0SMatthew G. Knepley PetscInt overlap2 = 1; 109c4762a1bSJed Brown 110*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 115c4762a1bSJed Brown 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeGeometryFVM(dm, &cell_geom, &face_geom)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(dm_view_geometry(dm, cell_geom, face_geom)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* redistribute */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm, &part)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-overlap2", &overlap2, NULL)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, overlap2, &redist_sf, &redist_dm)); 124c4762a1bSJed Brown if (redist_dm) { 1255f80ce2aSJacob Faibussowitsch CHKERRQ(redistribute_vec(redist_dm, redist_sf, &cell_geom)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(redistribute_vec(redist_dm, redist_sf, &face_geom)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) redist_sf, NULL, "-rd2_sf_view")); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "redistributed:\n")); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(dm_view_geometry(redist_dm, cell_geom, face_geom)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cell_geom)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&face_geom)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&redist_sf)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&redist_dm)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 137*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 138*b122ec5aSJacob Faibussowitsch return 0; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /*TEST 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: 0 145c4762a1bSJed Brown nsize: 3 146e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 8,1,1 -dm_plex_simplex 0 -dm_plex_adj_cone 1 -dm_plex_adj_closure 0 -petscpartitioner_type simple -dm_distribute_overlap 1 -overlap2 1 147c4762a1bSJed Brown 148c4762a1bSJed Brown TEST*/ 149