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