xref: /petsc/src/dm/impls/plex/tests/ex36.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
99371c9d4SSatish Balay static PetscErrorCode redistribute_vec(DM dist_dm, PetscSF sf, Vec *v) {
10c4762a1bSJed Brown   DM           dm, dist_v_dm;
11c4762a1bSJed Brown   PetscSection section, dist_section;
12c4762a1bSJed Brown   Vec          dist_v;
13c4762a1bSJed Brown   PetscMPIInt  rank, size, p;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(VecGetDM(*v, &dm));
179566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
189566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-rd_dm_view"));
199566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dist_dm, NULL, "-rd2_dm_view"));
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dist_v_dm));
229566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)*v), &dist_v));
239566063dSJacob Faibussowitsch   PetscCall(VecSetDM(dist_v, dist_v_dm));
249566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)*v), &dist_section));
259566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dist_v_dm, dist_section));
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)section, NULL, "-rd_section_view"));
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
30c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
31*48a46eb9SPierre Jolivet     if (p == rank) PetscCall(PetscObjectViewFromOptions((PetscObject)*v, NULL, "-rd_vec_view"));
329566063dSJacob Faibussowitsch     PetscCall(PetscBarrier((PetscObject)dm));
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
34c4762a1bSJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeField(dm, sf, section, *v, dist_section, dist_v));
36c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
37c4762a1bSJed Brown     if (p == rank) {
389566063dSJacob Faibussowitsch       PetscCall(PetscObjectViewFromOptions((PetscObject)dist_section, NULL, "-rd2_section_view"));
399566063dSJacob Faibussowitsch       PetscCall(PetscObjectViewFromOptions((PetscObject)dist_v, NULL, "-rd2_vec_view"));
40c4762a1bSJed Brown     }
419566063dSJacob Faibussowitsch     PetscCall(PetscBarrier((PetscObject)dm));
429566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&dist_section));
469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dist_v_dm));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(v));
49c4762a1bSJed Brown   *v = dist_v;
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
539371c9d4SSatish Balay static PetscErrorCode dm_view_geometry(DM dm, Vec cell_geom, Vec face_geom) {
54c4762a1bSJed Brown   DM                 cell_dm, face_dm;
55c4762a1bSJed Brown   PetscSection       cell_section, face_section;
56c4762a1bSJed Brown   const PetscScalar *cell_array, *face_array;
57c4762a1bSJed Brown   const PetscInt    *cells;
58c4762a1bSJed Brown   PetscInt           c, start_cell, end_cell;
59c4762a1bSJed Brown   PetscInt           f, start_face, end_face;
60c4762a1bSJed Brown   PetscInt           supportSize, offset;
61c4762a1bSJed Brown   PetscMPIInt        rank;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   PetscFunctionBegin;
649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* cells */
679566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &start_cell, &end_cell));
689566063dSJacob Faibussowitsch   PetscCall(VecGetDM(cell_geom, &cell_dm));
699566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cell_dm, &cell_section));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(cell_geom, &cell_array));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   for (c = start_cell; c < end_cell; ++c) {
73c4762a1bSJed Brown     const PetscFVCellGeom *geom;
749566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(cell_section, c, &offset));
75c4762a1bSJed Brown     geom = (PetscFVCellGeom *)&cell_array[offset];
7663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "rank %d c %" PetscInt_FMT " centroid %g,%g,%g vol %g\n", rank, c, (double)geom->centroid[0], (double)geom->centroid[1], (double)geom->centroid[2], (double)geom->volume));
77c4762a1bSJed Brown   }
789566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cell_geom, &cell_array));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* faces */
829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &start_face, &end_face));
839566063dSJacob Faibussowitsch   PetscCall(VecGetDM(face_geom, &face_dm));
849566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(face_dm, &face_section));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(face_geom, &face_array));
86c4762a1bSJed Brown   for (f = start_face; f < end_face; ++f) {
879566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, f, &cells));
889566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, f, &supportSize));
89c4762a1bSJed Brown     if (supportSize > 1) {
909566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(face_section, f, &offset));
9163a3b9bcSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "rank %d f %" PetscInt_FMT " cells %" PetscInt_FMT ",%" PetscInt_FMT " 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])));
92c4762a1bSJed Brown     }
93c4762a1bSJed Brown   }
949566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(cell_geom, &cell_array));
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
999371c9d4SSatish Balay int main(int argc, char **argv) {
10030602db0SMatthew G. Knepley   DM               dm, redist_dm;
101c4762a1bSJed Brown   PetscPartitioner part;
10230602db0SMatthew G. Knepley   PetscSF          redist_sf;
103c4762a1bSJed Brown   Vec              cell_geom, face_geom;
10430602db0SMatthew G. Knepley   PetscInt         overlap2 = 1;
105c4762a1bSJed Brown 
106327415f7SBarry Smith   PetscFunctionBeginUser;
1079566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1089566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
1099566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1119566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeGeometryFVM(dm, &cell_geom, &face_geom));
1149566063dSJacob Faibussowitsch   PetscCall(dm_view_geometry(dm, cell_geom, face_geom));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* redistribute */
1179566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm, &part));
1189566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part));
1199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-overlap2", &overlap2, NULL));
1209566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dm, overlap2, &redist_sf, &redist_dm));
121c4762a1bSJed Brown   if (redist_dm) {
1229566063dSJacob Faibussowitsch     PetscCall(redistribute_vec(redist_dm, redist_sf, &cell_geom));
1239566063dSJacob Faibussowitsch     PetscCall(redistribute_vec(redist_dm, redist_sf, &face_geom));
1249566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)redist_sf, NULL, "-rd2_sf_view"));
1259566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "redistributed:\n"));
1269566063dSJacob Faibussowitsch     PetscCall(dm_view_geometry(redist_dm, cell_geom, face_geom));
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cell_geom));
1309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&face_geom));
1319566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&redist_sf));
1329566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&redist_dm));
1339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
135b122ec5aSJacob Faibussowitsch   return 0;
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /*TEST
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   test:
141c4762a1bSJed Brown     suffix: 0
142c4762a1bSJed Brown     nsize: 3
143e600fa54SMatthew 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
144c4762a1bSJed Brown 
145c4762a1bSJed Brown TEST*/
146