xref: /petsc/src/dm/impls/plex/tests/ex70.c (revision 41e3e744d3c75c62dbad5fa29ebc2acfd1565676)
1*41e3e744Sksagiyam static char help[] = "Tests for DMPlexMarkBoundaryFaces()\n\n";
2*41e3e744Sksagiyam 
3*41e3e744Sksagiyam #include <petscdmplex.h>
4*41e3e744Sksagiyam #include <petscsf.h>
5*41e3e744Sksagiyam 
6*41e3e744Sksagiyam typedef struct {
7*41e3e744Sksagiyam   PetscInt overlap; /* The overlap size used when partitioning */
8*41e3e744Sksagiyam } AppCtx;
9*41e3e744Sksagiyam 
10*41e3e744Sksagiyam PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11*41e3e744Sksagiyam {
12*41e3e744Sksagiyam   PetscFunctionBegin;
13*41e3e744Sksagiyam   options->overlap = 0;
14*41e3e744Sksagiyam 
15*41e3e744Sksagiyam   PetscOptionsBegin(comm, "", "Options for DMPlexMarkBoundaryFaces() problem", "DMPLEX");
16*41e3e744Sksagiyam   PetscCall(PetscOptionsBoundedInt("-overlap", "The overlap size used", "ex70.c", options->overlap, &options->overlap, NULL, 0));
17*41e3e744Sksagiyam   PetscOptionsEnd();
18*41e3e744Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
19*41e3e744Sksagiyam }
20*41e3e744Sksagiyam 
21*41e3e744Sksagiyam static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
22*41e3e744Sksagiyam {
23*41e3e744Sksagiyam   PetscFunctionBeginUser;
24*41e3e744Sksagiyam   {
25*41e3e744Sksagiyam     const PetscInt faces[2] = {2, 2};
26*41e3e744Sksagiyam 
27*41e3e744Sksagiyam     PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, dm));
28*41e3e744Sksagiyam   }
29*41e3e744Sksagiyam   {
30*41e3e744Sksagiyam     PetscPartitioner part;
31*41e3e744Sksagiyam     PetscInt        *sizes  = NULL;
32*41e3e744Sksagiyam     PetscInt        *points = NULL;
33*41e3e744Sksagiyam     PetscMPIInt      rank, size;
34*41e3e744Sksagiyam 
35*41e3e744Sksagiyam     PetscCallMPI(MPI_Comm_rank(comm, &rank));
36*41e3e744Sksagiyam     PetscCallMPI(MPI_Comm_size(comm, &size));
37*41e3e744Sksagiyam     if (rank == 0) {
38*41e3e744Sksagiyam       PetscInt sizes1[2]  = {4, 4};
39*41e3e744Sksagiyam       PetscInt points1[8] = {3, 5, 6, 7, 0, 1, 2, 4};
40*41e3e744Sksagiyam 
41*41e3e744Sksagiyam       PetscCall(PetscMalloc2(2, &sizes, 8, &points));
42*41e3e744Sksagiyam       PetscCall(PetscArraycpy(sizes, sizes1, 2));
43*41e3e744Sksagiyam       PetscCall(PetscArraycpy(points, points1, 8));
44*41e3e744Sksagiyam     }
45*41e3e744Sksagiyam     PetscCall(DMPlexGetPartitioner(*dm, &part));
46*41e3e744Sksagiyam     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
47*41e3e744Sksagiyam     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
48*41e3e744Sksagiyam     PetscCall(PetscFree2(sizes, points));
49*41e3e744Sksagiyam   }
50*41e3e744Sksagiyam   {
51*41e3e744Sksagiyam     DM dmDist = NULL;
52*41e3e744Sksagiyam 
53*41e3e744Sksagiyam     PetscCall(DMSetAdjacency(*dm, -1, PETSC_FALSE, PETSC_TRUE));
54*41e3e744Sksagiyam     PetscCall(DMPlexDistribute(*dm, user->overlap, NULL, &dmDist));
55*41e3e744Sksagiyam     if (dmDist) {
56*41e3e744Sksagiyam       PetscCall(DMDestroy(dm));
57*41e3e744Sksagiyam       *dm = dmDist;
58*41e3e744Sksagiyam     }
59*41e3e744Sksagiyam   }
60*41e3e744Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
61*41e3e744Sksagiyam }
62*41e3e744Sksagiyam 
63*41e3e744Sksagiyam int main(int argc, char **argv)
64*41e3e744Sksagiyam {
65*41e3e744Sksagiyam   DM          dm;
66*41e3e744Sksagiyam   DMLabel     extLabel;
67*41e3e744Sksagiyam   MPI_Comm    comm;
68*41e3e744Sksagiyam   PetscMPIInt size;
69*41e3e744Sksagiyam   AppCtx      user;
70*41e3e744Sksagiyam 
71*41e3e744Sksagiyam   PetscFunctionBeginUser;
72*41e3e744Sksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
73*41e3e744Sksagiyam   comm = PETSC_COMM_WORLD;
74*41e3e744Sksagiyam   PetscCallMPI(MPI_Comm_size(comm, &size));
75*41e3e744Sksagiyam   if (size != 2) {
76*41e3e744Sksagiyam     PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 2.\n"));
77*41e3e744Sksagiyam     PetscCall(PetscFinalize());
78*41e3e744Sksagiyam     return 0;
79*41e3e744Sksagiyam   }
80*41e3e744Sksagiyam   PetscCall(ProcessOptions(comm, &user));
81*41e3e744Sksagiyam   PetscCall(CreateMesh(comm, &user, &dm));
82*41e3e744Sksagiyam   PetscCall(DMCreateLabel(dm, "exterior_facets"));
83*41e3e744Sksagiyam   PetscCall(DMGetLabel(dm, "exterior_facets", &extLabel));
84*41e3e744Sksagiyam   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, extLabel));
85*41e3e744Sksagiyam   PetscCall(PetscObjectSetName((PetscObject)dm, "Example_DM"));
86*41e3e744Sksagiyam   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
87*41e3e744Sksagiyam   PetscCall(DMDestroy(&dm));
88*41e3e744Sksagiyam   PetscCall(PetscFinalize());
89*41e3e744Sksagiyam   return 0;
90*41e3e744Sksagiyam }
91*41e3e744Sksagiyam 
92*41e3e744Sksagiyam /*TEST
93*41e3e744Sksagiyam 
94*41e3e744Sksagiyam   test:
95*41e3e744Sksagiyam     nsize: 2
96*41e3e744Sksagiyam     requires: triangle
97*41e3e744Sksagiyam     args: -overlap {{0 1}separate output} -dm_view ascii::ascii_info_detail
98*41e3e744Sksagiyam 
99*41e3e744Sksagiyam TEST*/
100