xref: /petsc/src/dm/impls/plex/tests/ex29.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 static char help[] = "Test scalable partitioning on distributed meshes\n\n";
2 
3 #include <petscdmplex.h>
4 
5 enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_OVERLAP};
6 
7 typedef struct {
8   PetscLogEvent createMeshEvent;
9   PetscLogStage stages[4];
10   /* Domain and mesh definition */
11   PetscInt overlap; /* The cell overlap to use during partitioning */
12 } AppCtx;
13 
14 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15 {
16   PetscFunctionBegin;
17   options->overlap = PETSC_FALSE;
18 
19   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
20   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL,0));
21   PetscOptionsEnd();
22 
23   PetscCall(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent));
24   PetscCall(PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]));
25   PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
26   PetscCall(PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]));
27   PetscCall(PetscLogStageRegister("MeshOverlap",    &options->stages[STAGE_OVERLAP]));
28   PetscFunctionReturn(0);
29 }
30 
31 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
32 {
33   PetscMPIInt    rank, size;
34 
35   PetscFunctionBegin;
36   PetscCall(PetscLogEventBegin(user->createMeshEvent,0,0,0,0));
37   PetscCallMPI(MPI_Comm_rank(comm, &rank));
38   PetscCallMPI(MPI_Comm_size(comm, &size));
39   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
40   PetscCall(DMCreate(comm, dm));
41   PetscCall(DMSetType(*dm, DMPLEX));
42   PetscCall(DMSetFromOptions(*dm));
43   PetscCall(PetscLogStagePop());
44   {
45     DM               pdm = NULL;
46     PetscPartitioner part;
47 
48     PetscCall(DMPlexGetPartitioner(*dm, &part));
49     PetscCall(PetscPartitionerSetFromOptions(part));
50     /* Distribute mesh over processes */
51     PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
52     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
53     if (pdm) {
54       PetscCall(DMDestroy(dm));
55       *dm  = pdm;
56     }
57     PetscCall(PetscLogStagePop());
58   }
59   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
60   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_"));
61   PetscCall(DMSetFromOptions(*dm));
62   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, ""));
63   PetscCall(PetscLogStagePop());
64   if (user->overlap) {
65     DM odm = NULL;
66     /* Add the level-1 overlap to refined mesh */
67     PetscCall(PetscLogStagePush(user->stages[STAGE_OVERLAP]));
68     PetscCall(DMPlexDistributeOverlap(*dm, 1, NULL, &odm));
69     if (odm) {
70       PetscCall(DMView(odm, PETSC_VIEWER_STDOUT_WORLD));
71       PetscCall(DMDestroy(dm));
72       *dm = odm;
73     }
74     PetscCall(PetscLogStagePop());
75   }
76   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
77   PetscCall(PetscLogEventEnd(user->createMeshEvent,0,0,0,0));
78   PetscFunctionReturn(0);
79 }
80 
81 int main(int argc, char **argv)
82 {
83   DM             dm, pdm;
84   AppCtx         user;                 /* user-defined work context */
85   PetscPartitioner part;
86 
87   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
88   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
89   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
90   PetscCall(DMPlexGetPartitioner(dm, &part));
91   PetscCall(PetscPartitionerSetFromOptions(part));
92   PetscCall(DMPlexDistribute(dm, user.overlap, NULL, &pdm));
93   if (pdm) PetscCall(DMViewFromOptions(pdm, NULL, "-pdm_view"));
94   PetscCall(DMDestroy(&dm));
95   PetscCall(DMDestroy(&pdm));
96   PetscCall(PetscFinalize());
97   return 0;
98 }
99 
100 /*TEST
101 
102   test:
103     suffix: 0
104     requires: ctetgen
105     args: -dm_plex_dim 3 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
106   test:
107     suffix: 1
108     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
109   test:
110     suffix: quad_0
111     nsize: 2
112     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view
113 
114 TEST*/
115