xref: /petsc/src/dm/impls/plex/tests/ex29.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscFunctionBegin;
17c4762a1bSJed Brown   options->overlap = PETSC_FALSE;
18c4762a1bSJed Brown 
19d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL,0));
21d0609cedSBarry Smith   PetscOptionsEnd();
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent));
249566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]));
259566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
269566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]));
279566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MeshOverlap",    &options->stages[STAGE_OVERLAP]));
28c4762a1bSJed Brown   PetscFunctionReturn(0);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
31c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
32c4762a1bSJed Brown {
33c4762a1bSJed Brown   PetscMPIInt    rank, size;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(user->createMeshEvent,0,0,0,0));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
399566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD]));
409566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
419566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
429566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
439566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
44c4762a1bSJed Brown   {
45c4762a1bSJed Brown     DM               pdm = NULL;
46c4762a1bSJed Brown     PetscPartitioner part;
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
499566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
50c4762a1bSJed Brown     /* Distribute mesh over processes */
519566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
529566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
53c4762a1bSJed Brown     if (pdm) {
549566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
55c4762a1bSJed Brown       *dm  = pdm;
56c4762a1bSJed Brown     }
579566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
58c4762a1bSJed Brown   }
599566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE]));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_"));
619566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, ""));
639566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
64c4762a1bSJed Brown   if (user->overlap) {
65c4762a1bSJed Brown     DM odm = NULL;
66c4762a1bSJed Brown     /* Add the level-1 overlap to refined mesh */
679566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(user->stages[STAGE_OVERLAP]));
689566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeOverlap(*dm, 1, NULL, &odm));
69c4762a1bSJed Brown     if (odm) {
709566063dSJacob Faibussowitsch       PetscCall(DMView(odm, PETSC_VIEWER_STDOUT_WORLD));
719566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
72c4762a1bSJed Brown       *dm = odm;
73c4762a1bSJed Brown     }
749566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
75c4762a1bSJed Brown   }
769566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(user->createMeshEvent,0,0,0,0));
78c4762a1bSJed Brown   PetscFunctionReturn(0);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown int main(int argc, char **argv)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   DM             dm, pdm;
84c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
85c4762a1bSJed Brown   PetscPartitioner part;
86c4762a1bSJed Brown 
87*327415f7SBarry Smith   PetscFunctionBeginUser;
889566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
899566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
909566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm, &part));
929566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part));
939566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dm, user.overlap, NULL, &pdm));
949566063dSJacob Faibussowitsch   if (pdm) PetscCall(DMViewFromOptions(pdm, NULL, "-pdm_view"));
959566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
969566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pdm));
979566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
98b122ec5aSJacob Faibussowitsch   return 0;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*TEST
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   test:
104c4762a1bSJed Brown     suffix: 0
105c4762a1bSJed Brown     requires: ctetgen
10630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
107c4762a1bSJed Brown   test:
108c4762a1bSJed Brown     suffix: 1
10930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
110c4762a1bSJed Brown   test:
111c4762a1bSJed Brown     suffix: quad_0
112c4762a1bSJed Brown     nsize: 2
11330602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view
114c4762a1bSJed Brown 
115c4762a1bSJed Brown TEST*/
116