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; 17*9566063dSJacob Faibussowitsch PetscCall(VecGetDM(*v, &dm)); 18*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 19*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-rd_dm_view")); 20*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dist_dm, NULL, "-rd2_dm_view")); 21c4762a1bSJed Brown 22*9566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dist_v_dm)); 23*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject) *v), &dist_v)); 24*9566063dSJacob Faibussowitsch PetscCall(VecSetDM(dist_v, dist_v_dm)); 25*9566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) *v), &dist_section)); 26*9566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dist_v_dm, dist_section)); 27c4762a1bSJed Brown 28*9566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) section, NULL, "-rd_section_view")); 29*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank)); 30*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 31c4762a1bSJed Brown for (p = 0; p < size; ++p) { 32c4762a1bSJed Brown if (p == rank) { 33*9566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) *v, NULL, "-rd_vec_view"));} 34*9566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject) dm)); 35*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 36c4762a1bSJed Brown } 37*9566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeField(dm, sf, section, *v, dist_section, dist_v)); 38c4762a1bSJed Brown for (p = 0; p < size; ++p) { 39c4762a1bSJed Brown if (p == rank) { 40*9566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) dist_section, NULL, "-rd2_section_view")); 41*9566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) dist_v, NULL, "-rd2_vec_view")); 42c4762a1bSJed Brown } 43*9566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject) dm)); 44*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47*9566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&dist_section)); 48*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dist_v_dm)); 49c4762a1bSJed Brown 50*9566063dSJacob Faibussowitsch PetscCall(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; 67*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* cells */ 70*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &start_cell, &end_cell)); 71*9566063dSJacob Faibussowitsch PetscCall(VecGetDM(cell_geom, &cell_dm)); 72*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cell_dm, &cell_section)); 73*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cell_geom, &cell_array)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown for (c = start_cell; c < end_cell; ++c) { 76c4762a1bSJed Brown const PetscFVCellGeom *geom; 77*9566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(cell_section, c, &offset)); 78c4762a1bSJed Brown geom = (PetscFVCellGeom*)&cell_array[offset]; 79*9566063dSJacob Faibussowitsch PetscCall(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 } 81*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 82*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cell_geom, &cell_array)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* faces */ 85*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &start_face, &end_face)); 86*9566063dSJacob Faibussowitsch PetscCall(VecGetDM(face_geom, &face_dm)); 87*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(face_dm, &face_section)); 88*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(face_geom, &face_array)); 89c4762a1bSJed Brown for (f = start_face; f < end_face; ++f) { 90*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, f, &cells)); 91*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, f, &supportSize)); 92c4762a1bSJed Brown if (supportSize > 1) { 93*9566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(face_section, f, &offset)); 94*9566063dSJacob Faibussowitsch PetscCall(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 } 97*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 98*9566063dSJacob Faibussowitsch PetscCall(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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 111*9566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 112*9566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 113*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 114*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 115c4762a1bSJed Brown 116*9566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGeometryFVM(dm, &cell_geom, &face_geom)); 117*9566063dSJacob Faibussowitsch PetscCall(dm_view_geometry(dm, cell_geom, face_geom)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* redistribute */ 120*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 121*9566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 122*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-overlap2", &overlap2, NULL)); 123*9566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, overlap2, &redist_sf, &redist_dm)); 124c4762a1bSJed Brown if (redist_dm) { 125*9566063dSJacob Faibussowitsch PetscCall(redistribute_vec(redist_dm, redist_sf, &cell_geom)); 126*9566063dSJacob Faibussowitsch PetscCall(redistribute_vec(redist_dm, redist_sf, &face_geom)); 127*9566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) redist_sf, NULL, "-rd2_sf_view")); 128*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "redistributed:\n")); 129*9566063dSJacob Faibussowitsch PetscCall(dm_view_geometry(redist_dm, cell_geom, face_geom)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cell_geom)); 133*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&face_geom)); 134*9566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&redist_sf)); 135*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&redist_dm)); 136*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 137*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 138b122ec5aSJacob 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