xref: /petsc/src/dm/impls/plex/tests/ex29.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 static char help[] = "Test scalalbe 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  dim;                             /* The topological mesh dimension */
12   PetscInt  faces[3];                        /* Number of faces per dimension */
13   PetscBool simplex;                         /* Use simplices or hexes */
14   char      filename[PETSC_MAX_PATH_LEN];    /* Import mesh from file */
15   PetscInt  overlap;                         /* The cell overlap to use during partitioning */
16 } AppCtx;
17 
18 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
19 {
20   PetscInt        dim;
21   PetscErrorCode  ierr;
22 
23   PetscFunctionBegin;
24   options->dim         = 2;
25   options->simplex     = PETSC_TRUE;
26   options->filename[0] = '\0';
27   options->overlap     = PETSC_FALSE;
28   options->faces[0]    = 1;
29   options->faces[1]    = 1;
30   options->faces[2]    = 1;
31 
32   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
33   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex27.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
34   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex27.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
35   ierr = PetscOptionsString("-filename", "The mesh file", "", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
36   ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex27.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr);
37   dim = options->dim;
38   ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex27.c", options->faces, &dim, NULL);CHKERRQ(ierr);
39   if (dim) options->dim = dim;
40   ierr = PetscOptionsEnd();
41 
42   ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr);
43   ierr = PetscLogStageRegister("MeshLoad",       &options->stages[STAGE_LOAD]);CHKERRQ(ierr);
44   ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr);
45   ierr = PetscLogStageRegister("MeshRefine",     &options->stages[STAGE_REFINE]);CHKERRQ(ierr);
46   ierr = PetscLogStageRegister("MeshOverlap",    &options->stages[STAGE_OVERLAP]);CHKERRQ(ierr);
47   PetscFunctionReturn(0);
48 }
49 
50 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
51 {
52   PetscInt       dim      = user->dim;
53   PetscInt      *faces    = user->faces;
54   PetscBool      simplex  = user->simplex;
55   const char    *filename = user->filename;
56   size_t         len;
57   PetscMPIInt    rank, size;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
62   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
63   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
64   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
65   ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr);
66   if (len) {
67     ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);
68   } else {
69     ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
70   }
71   ierr = PetscLogStagePop();CHKERRQ(ierr);
72   {
73     DM               pdm = NULL;
74     PetscPartitioner part;
75 
76     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
77     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
78     /* Distribute mesh over processes */
79     ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr);
80     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
81     if (pdm) {
82       ierr = DMDestroy(dm);CHKERRQ(ierr);
83       *dm  = pdm;
84     }
85     ierr = PetscLogStagePop();CHKERRQ(ierr);
86   }
87   ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr);
88   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
89   ierr = PetscLogStagePop();CHKERRQ(ierr);
90   if (user->overlap) {
91     DM odm = NULL;
92     /* Add the level-1 overlap to refined mesh */
93     ierr = PetscLogStagePush(user->stages[STAGE_OVERLAP]);CHKERRQ(ierr);
94     ierr = DMPlexDistributeOverlap(*dm, 1, NULL, &odm);CHKERRQ(ierr);
95     if (odm) {
96       ierr = DMView(odm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
97       ierr = DMDestroy(dm);CHKERRQ(ierr);
98       *dm = odm;
99     }
100     ierr = PetscLogStagePop();CHKERRQ(ierr);
101   }
102   if (simplex){ierr = PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh");CHKERRQ(ierr);}
103   else{ierr = PetscObjectSetName((PetscObject) *dm, "Tensor Product Mesh");CHKERRQ(ierr);}
104   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
105   ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 int main(int argc, char **argv)
110 {
111   DM             dm, pdm;
112   AppCtx         user;                 /* user-defined work context */
113   PetscErrorCode ierr;
114   PetscPartitioner part;
115 
116   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
117   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
118   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
119   ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr);
120   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
121   ierr = DMPlexDistribute(dm, user.overlap, NULL, &pdm);CHKERRQ(ierr);
122   if (pdm) {ierr = DMViewFromOptions(pdm, NULL, "-pdm_view");CHKERRQ(ierr);}
123   ierr = DMDestroy(&dm);CHKERRQ(ierr);
124   ierr = DMDestroy(&pdm);CHKERRQ(ierr);
125   ierr = PetscFinalize();
126   return ierr;
127 }
128 
129 /*TEST
130 
131   test:
132     suffix: 0
133     requires: ctetgen
134     args: -dim 3 -dm_refine 2 -petscpartitioner_type simple -dm_view
135   test:
136     suffix: 1
137     args: -dim 3 -simplex 0 -dm_refine 2 -petscpartitioner_type simple -dm_view
138   test:
139     suffix: quad_0
140     nsize: 2
141     args: -dim 3 -simplex 0 -dm_refine 2 -petscpartitioner_type simple -dm_view -pdm_view
142 
143 TEST*/
144