xref: /petsc/src/dm/impls/plex/tests/ex29.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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);
225f80ce2aSJacob 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 
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]));
295f80ce2aSJacob 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;
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(user->createMeshEvent,0,0,0,0));
395f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
405f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_LOAD]));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
46c4762a1bSJed Brown   {
47c4762a1bSJed Brown     DM               pdm = NULL;
48c4762a1bSJed Brown     PetscPartitioner part;
49c4762a1bSJed Brown 
505f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(*dm, &part));
515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetFromOptions(part));
52c4762a1bSJed Brown     /* Distribute mesh over processes */
535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]));
545f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistribute(*dm, 0, NULL, &pdm));
55c4762a1bSJed Brown     if (pdm) {
565f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(dm));
57c4762a1bSJed Brown       *dm  = pdm;
58c4762a1bSJed Brown     }
595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
60c4762a1bSJed Brown   }
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(user->stages[STAGE_REFINE]));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_"));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, ""));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
66c4762a1bSJed Brown   if (user->overlap) {
67c4762a1bSJed Brown     DM odm = NULL;
68c4762a1bSJed Brown     /* Add the level-1 overlap to refined mesh */
695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(user->stages[STAGE_OVERLAP]));
705f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistributeOverlap(*dm, 1, NULL, &odm));
71c4762a1bSJed Brown     if (odm) {
725f80ce2aSJacob Faibussowitsch       CHKERRQ(DMView(odm, PETSC_VIEWER_STDOUT_WORLD));
735f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(dm));
74c4762a1bSJed Brown       *dm = odm;
75c4762a1bSJed Brown     }
765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
77c4762a1bSJed Brown   }
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
795f80ce2aSJacob 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   PetscPartitioner part;
88c4762a1bSJed Brown 
89*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetPartitioner(dm, &part));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerSetFromOptions(part));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(dm, user.overlap, NULL, &pdm));
955f80ce2aSJacob Faibussowitsch   if (pdm) CHKERRQ(DMViewFromOptions(pdm, NULL, "-pdm_view"));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&pdm));
98*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
99*b122ec5aSJacob 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