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