xref: /petsc/src/dm/impls/plex/tests/ex29.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
21c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex29.c", options->overlap, &options->overlap, NULL,0));
231e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
24c4762a1bSJed Brown 
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(user->createMeshEvent,0,0,0,0));
39*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
40*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_LOAD]));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
46c4762a1bSJed Brown   {
47c4762a1bSJed Brown     DM               pdm = NULL;
48c4762a1bSJed Brown     PetscPartitioner part;
49c4762a1bSJed Brown 
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(*dm, &part));
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetFromOptions(part));
52c4762a1bSJed Brown     /* Distribute mesh over processes */
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistribute(*dm, 0, NULL, &pdm));
55c4762a1bSJed Brown     if (pdm) {
56*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(dm));
57c4762a1bSJed Brown       *dm  = pdm;
58c4762a1bSJed Brown     }
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
60c4762a1bSJed Brown   }
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_REFINE]));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_"));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, ""));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
66c4762a1bSJed Brown   if (user->overlap) {
67c4762a1bSJed Brown     DM odm = NULL;
68c4762a1bSJed Brown     /* Add the level-1 overlap to refined mesh */
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(user->stages[STAGE_OVERLAP]));
70*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistributeOverlap(*dm, 1, NULL, &odm));
71c4762a1bSJed Brown     if (odm) {
72*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMView(odm, PETSC_VIEWER_STDOUT_WORLD));
73*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(dm));
74c4762a1bSJed Brown       *dm = odm;
75c4762a1bSJed Brown     }
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
77c4762a1bSJed Brown   }
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode ierr;
88c4762a1bSJed Brown   PetscPartitioner part;
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm, &part));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm, user.overlap, NULL, &pdm));
96*5f80ce2aSJacob Faibussowitsch   if (pdm) CHKERRQ(DMViewFromOptions(pdm, NULL, "-pdm_view"));
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&pdm));
99c4762a1bSJed Brown   ierr = PetscFinalize();
100c4762a1bSJed Brown   return ierr;
101c4762a1bSJed Brown }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown /*TEST
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   test:
106c4762a1bSJed Brown     suffix: 0
107c4762a1bSJed Brown     requires: ctetgen
10830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
109c4762a1bSJed Brown   test:
110c4762a1bSJed Brown     suffix: 1
11130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view
112c4762a1bSJed Brown   test:
113c4762a1bSJed Brown     suffix: quad_0
114c4762a1bSJed Brown     nsize: 2
11530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -post_dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view
116c4762a1bSJed Brown 
117c4762a1bSJed Brown TEST*/
118