1f5a72eb3SToby Isaac const char help[] = "Set up a PetscSF for halo exchange between local vectors";
2f5a72eb3SToby Isaac
3f5a72eb3SToby Isaac #include <petscdmplex.h>
4f5a72eb3SToby Isaac #include <petscsf.h>
5f5a72eb3SToby Isaac
main(int argc,char ** argv)6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7*d71ae5a4SJacob Faibussowitsch {
8f5a72eb3SToby Isaac DM dm;
9f5a72eb3SToby Isaac Vec u;
10f5a72eb3SToby Isaac PetscSection section;
11f5a72eb3SToby Isaac PetscInt dim, numFields, numBC, i;
12f5a72eb3SToby Isaac PetscMPIInt rank;
13f5a72eb3SToby Isaac PetscInt numComp[2];
14f5a72eb3SToby Isaac PetscInt numDof[12];
15f5a72eb3SToby Isaac PetscInt *remoteOffsets;
16f5a72eb3SToby Isaac PetscSF pointSF;
17f5a72eb3SToby Isaac PetscSF sectionSF;
18f5a72eb3SToby Isaac PetscScalar *array;
19f5a72eb3SToby Isaac
20327415f7SBarry Smith PetscFunctionBeginUser;
21f5a72eb3SToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help));
22f5a72eb3SToby Isaac /* Create a mesh */
23f5a72eb3SToby Isaac PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
24f5a72eb3SToby Isaac PetscCall(DMSetType(dm, DMPLEX));
25f5a72eb3SToby Isaac PetscCall(DMSetFromOptions(dm));
26f5a72eb3SToby Isaac PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
27f5a72eb3SToby Isaac PetscCall(DMGetDimension(dm, &dim));
28f5a72eb3SToby Isaac
29f5a72eb3SToby Isaac /** Describe the solution variables that are discretized on the mesh */
30f5a72eb3SToby Isaac /* Create scalar field u and a vector field v */
31f5a72eb3SToby Isaac numFields = 2;
32f5a72eb3SToby Isaac numComp[0] = 1;
33f5a72eb3SToby Isaac numComp[1] = dim;
34f5a72eb3SToby Isaac for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0;
35f5a72eb3SToby Isaac /* Let u be defined on cells */
36f5a72eb3SToby Isaac numDof[0 * (dim + 1) + dim] = 1;
37f5a72eb3SToby Isaac /* Let v be defined on vertices */
38f5a72eb3SToby Isaac numDof[1 * (dim + 1) + 0] = dim;
39f5a72eb3SToby Isaac /* No boundary conditions */
40f5a72eb3SToby Isaac numBC = 0;
41f5a72eb3SToby Isaac
42f5a72eb3SToby Isaac /** Create a PetscSection to handle the layout of the discretized variables */
43f5a72eb3SToby Isaac PetscCall(DMSetNumFields(dm, numFields));
44f5a72eb3SToby Isaac PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, §ion));
45f5a72eb3SToby Isaac /* Name the Field variables */
46f5a72eb3SToby Isaac PetscCall(PetscSectionSetFieldName(section, 0, "u"));
47f5a72eb3SToby Isaac PetscCall(PetscSectionSetFieldName(section, 1, "v"));
48f5a72eb3SToby Isaac /* Tell the DM to use this data layout */
49f5a72eb3SToby Isaac PetscCall(DMSetLocalSection(dm, section));
50f5a72eb3SToby Isaac
51f5a72eb3SToby Isaac /** Construct the communication pattern for halo exchange between local vectors */
52f5a72eb3SToby Isaac /* Get the point SF: an object that says which copies of mesh points (cells,
53f5a72eb3SToby Isaac * vertices, faces, edges) are copies of points on other processes */
54f5a72eb3SToby Isaac PetscCall(DMGetPointSF(dm, &pointSF));
55f5a72eb3SToby Isaac /* Relate the locations of ghost degrees of freedom on this process
56f5a72eb3SToby Isaac * to their locations of the non-ghost copies on a different process */
57f5a72eb3SToby Isaac PetscCall(PetscSFCreateRemoteOffsets(pointSF, section, section, &remoteOffsets));
58f5a72eb3SToby Isaac /* Use that information to construct a star forest for halo exchange
59f5a72eb3SToby Isaac * for data described by the local section */
60f5a72eb3SToby Isaac PetscCall(PetscSFCreateSectionSF(pointSF, section, remoteOffsets, section, §ionSF));
61f5a72eb3SToby Isaac PetscCall(PetscFree(remoteOffsets));
62f5a72eb3SToby Isaac
63f5a72eb3SToby Isaac /** Demo of halo exchange */
64f5a72eb3SToby Isaac /* Create a Vec with this layout */
65f5a72eb3SToby Isaac PetscCall(DMCreateLocalVector(dm, &u));
66f5a72eb3SToby Isaac PetscCall(PetscObjectSetName((PetscObject)u, "Local vector"));
67f5a72eb3SToby Isaac PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
68f5a72eb3SToby Isaac /* Set all mesh values to the MPI rank */
69f5a72eb3SToby Isaac PetscCall(VecSet(u, (PetscScalar)rank));
70f5a72eb3SToby Isaac /* Get the raw array of values */
71f5a72eb3SToby Isaac PetscCall(VecGetArrayWrite(u, &array));
72f5a72eb3SToby Isaac /*** HALO EXCHANGE ***/
73f5a72eb3SToby Isaac PetscCall(PetscSFBcastBegin(sectionSF, MPIU_SCALAR, array, array, MPI_REPLACE));
74f5a72eb3SToby Isaac /* local work can be done between Begin() and End() */
75f5a72eb3SToby Isaac PetscCall(PetscSFBcastEnd(sectionSF, MPIU_SCALAR, array, array, MPI_REPLACE));
76f5a72eb3SToby Isaac /* Restore the raw array of values */
77f5a72eb3SToby Isaac PetscCall(VecRestoreArrayWrite(u, &array));
78f5a72eb3SToby Isaac /* View the results: should show which process has the non-ghost copy of each degree of freedom */
79f5a72eb3SToby Isaac PetscCall(PetscSectionVecView(section, u, PETSC_VIEWER_STDOUT_WORLD));
80f5a72eb3SToby Isaac PetscCall(VecDestroy(&u));
81f5a72eb3SToby Isaac
82f5a72eb3SToby Isaac /** Cleanup */
83f5a72eb3SToby Isaac PetscCall(PetscSFDestroy(§ionSF));
84f5a72eb3SToby Isaac PetscCall(PetscSectionDestroy(§ion));
85f5a72eb3SToby Isaac PetscCall(DMDestroy(&dm));
86f5a72eb3SToby Isaac PetscCall(PetscFinalize());
87f5a72eb3SToby Isaac return 0;
88f5a72eb3SToby Isaac }
89f5a72eb3SToby Isaac
90f5a72eb3SToby Isaac /*TEST
91f5a72eb3SToby Isaac
92f5a72eb3SToby Isaac # Test on a 1D mesh with overlap
93f5a72eb3SToby Isaac test:
94f5a72eb3SToby Isaac nsize: 3
95f5a72eb3SToby Isaac requires: !complex
96f5a72eb3SToby Isaac args: -dm_plex_dim 1 -dm_plex_box_faces 3 -dm_refine_pre 1 -petscpartitioner_type simple -dm_distribute_overlap 1
97f5a72eb3SToby Isaac
98f5a72eb3SToby Isaac TEST*/
99