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