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 PetscErrorCode ierr; 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscFunctionBegin; 18c4762a1bSJed Brown ierr = VecGetDM(*v, &dm);CHKERRQ(ierr); 19c4762a1bSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 20c4762a1bSJed Brown ierr = DMViewFromOptions(dm, NULL, "-rd_dm_view");CHKERRQ(ierr); 21c4762a1bSJed Brown ierr = DMViewFromOptions(dist_dm, NULL, "-rd2_dm_view");CHKERRQ(ierr); 22c4762a1bSJed Brown 23c4762a1bSJed Brown ierr = DMClone(dm, &dist_v_dm);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = VecCreate(PetscObjectComm((PetscObject) *v), &dist_v);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = VecSetDM(dist_v, dist_v_dm);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscSectionCreate(PetscObjectComm((PetscObject) *v), &dist_section);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = DMSetLocalSection(dist_v_dm, dist_section);CHKERRQ(ierr); 28c4762a1bSJed Brown 29c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) section, NULL, "-rd_section_view");CHKERRQ(ierr); 30*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRMPI(ierr); 31*ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 32c4762a1bSJed Brown for (p = 0; p < size; ++p) { 33c4762a1bSJed Brown if (p == rank) { 34c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) *v, NULL, "-rd_vec_view");CHKERRQ(ierr);} 35c4762a1bSJed Brown ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown ierr = DMPlexDistributeField(dm, sf, section, *v, dist_section, dist_v);CHKERRQ(ierr); 39c4762a1bSJed Brown for (p = 0; p < size; ++p) { 40c4762a1bSJed Brown if (p == rank) { 41c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) dist_section, NULL, "-rd2_section_view");CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) dist_v, NULL, "-rd2_vec_view");CHKERRQ(ierr); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown ierr = PetscSectionDestroy(&dist_section);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = DMDestroy(&dist_v_dm);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = VecDestroy(v);CHKERRQ(ierr); 52c4762a1bSJed Brown *v = dist_v; 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown static PetscErrorCode dm_view_geometry(DM dm, Vec cell_geom, Vec face_geom) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown DM cell_dm, face_dm; 59c4762a1bSJed Brown PetscSection cell_section, face_section; 60c4762a1bSJed Brown const PetscScalar *cell_array, *face_array; 61c4762a1bSJed Brown const PetscInt *cells; 62c4762a1bSJed Brown PetscInt c, start_cell, end_cell; 63c4762a1bSJed Brown PetscInt f, start_face, end_face; 64c4762a1bSJed Brown PetscInt supportSize, offset; 65c4762a1bSJed Brown PetscMPIInt rank; 66c4762a1bSJed Brown PetscErrorCode ierr; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBegin; 69*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* cells */ 72c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &start_cell, &end_cell);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = VecGetDM(cell_geom, &cell_dm);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = DMGetLocalSection(cell_dm, &cell_section);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = VecGetArrayRead(cell_geom, &cell_array);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown for (c = start_cell; c < end_cell; ++c) { 78c4762a1bSJed Brown const PetscFVCellGeom *geom; 79c4762a1bSJed Brown ierr = PetscSectionGetOffset(cell_section, c, &offset);CHKERRQ(ierr); 80c4762a1bSJed Brown geom = (PetscFVCellGeom*)&cell_array[offset]; 81c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = VecRestoreArrayRead(cell_geom, &cell_array);CHKERRQ(ierr); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* faces */ 87c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 1, &start_face, &end_face);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = VecGetDM(face_geom, &face_dm);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = DMGetLocalSection(face_dm, &face_section);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = VecGetArrayRead(face_geom, &face_array);CHKERRQ(ierr); 91c4762a1bSJed Brown for (f = start_face; f < end_face; ++f) { 92c4762a1bSJed Brown ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 94c4762a1bSJed Brown if (supportSize > 1) { 95c4762a1bSJed Brown ierr = PetscSectionGetOffset(face_section, f, &offset);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = 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]));CHKERRQ(ierr); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown } 99c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = VecRestoreArrayRead(cell_geom, &cell_array);CHKERRQ(ierr); 101c4762a1bSJed Brown PetscFunctionReturn(0); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown int main(int argc, char **argv) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown DM dm, dist_dm, redist_dm; 107c4762a1bSJed Brown PetscPartitioner part; 108c4762a1bSJed Brown PetscSF dist_sf, redist_sf; 109c4762a1bSJed Brown Vec cell_geom, face_geom; 110c4762a1bSJed Brown PetscInt overlap = 1, overlap2 = 1; 111c4762a1bSJed Brown PetscMPIInt rank; 112c4762a1bSJed Brown const char *filename = "gminc_1d.exo"; 113c4762a1bSJed Brown PetscErrorCode ierr; 114c4762a1bSJed Brown 115c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 116*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 117c4762a1bSJed Brown if (0) { 118c4762a1bSJed Brown ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);CHKERRQ(ierr); 119c4762a1bSJed Brown } else { 120c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 3, PETSC_FALSE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown ierr = DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL, NULL, "-overlap", &overlap, NULL);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = DMPlexDistribute(dm, overlap, &dist_sf, &dist_dm);CHKERRQ(ierr); 128c4762a1bSJed Brown if (dist_dm) { 129c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 130c4762a1bSJed Brown dm = dist_dm; 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown ierr = DMPlexComputeGeometryFVM(dm, &cell_geom, &face_geom);CHKERRQ(ierr); 134c4762a1bSJed Brown 135c4762a1bSJed Brown ierr = dm_view_geometry(dm, cell_geom, face_geom);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* redistribute */ 138c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL, NULL, "-overlap2", &overlap2, NULL);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMPlexDistribute(dm, overlap2, &redist_sf, &redist_dm);CHKERRQ(ierr); 142c4762a1bSJed Brown if (redist_dm) { 143c4762a1bSJed Brown ierr = redistribute_vec(redist_dm, redist_sf, &cell_geom);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = redistribute_vec(redist_dm, redist_sf, &face_geom);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) redist_sf, NULL, "-rd2_sf_view");CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "redistributed:\n");CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = dm_view_geometry(redist_dm, cell_geom, face_geom);CHKERRQ(ierr); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown ierr = VecDestroy(&cell_geom);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = VecDestroy(&face_geom);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscSFDestroy(&dist_sf);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = PetscSFDestroy(&redist_sf);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = DMDestroy(&redist_dm);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = PetscFinalize(); 157c4762a1bSJed Brown return ierr; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /*TEST 161c4762a1bSJed Brown 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown suffix: 0 164c4762a1bSJed Brown nsize: 3 165c4762a1bSJed Brown args: -overlap 1 -overlap2 1 -dm_plex_box_faces 8,1,1 -petscpartitioner_type simple 166c4762a1bSJed Brown 167c4762a1bSJed Brown TEST*/ 168