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