xref: /petsc/src/dm/impls/plex/tests/ex13.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Orient a mesh in parallel\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   /* Domain and mesh definition */
7c4762a1bSJed Brown   PetscBool testPartition; /* Use a fixed partitioning for testing */
8c4762a1bSJed Brown   PetscInt  testNum;       /* Labels the different test partitions */
9c4762a1bSJed Brown } AppCtx;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscFunctionBeginUser;
14c4762a1bSJed Brown   options->testPartition = PETSC_TRUE;
15c4762a1bSJed Brown   options->testNum       = 0;
16c4762a1bSJed Brown 
17d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL,0));
20d0609cedSBarry Smith   PetscOptionsEnd();
21c4762a1bSJed Brown   PetscFunctionReturn(0);
22c4762a1bSJed Brown }
23c4762a1bSJed Brown 
24c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   DM             dmDist = NULL;
2730602db0SMatthew G. Knepley   PetscBool      simplex;
2830602db0SMatthew G. Knepley   PetscInt       dim;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   PetscFunctionBeginUser;
319566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
329566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
339566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
349566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
359566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
369566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(*dm, &simplex));
37c4762a1bSJed Brown   if (user->testPartition) {
38c4762a1bSJed Brown     PetscPartitioner part;
39c4762a1bSJed Brown     PetscInt        *sizes  = NULL;
40c4762a1bSJed Brown     PetscInt        *points = NULL;
41c4762a1bSJed Brown     PetscMPIInt      rank, size;
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
449566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
45dd400576SPatrick Sanan     if (rank == 0) {
4630602db0SMatthew G. Knepley       if (dim == 2 && simplex && size == 2) {
47c4762a1bSJed Brown         switch (user->testNum) {
48c4762a1bSJed Brown         case 0: {
49c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {4, 4};
50c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4};
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 8, &points));
539566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
549566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, triPoints_p2, 8));
55c4762a1bSJed Brown           break;}
56c4762a1bSJed Brown         case 1: {
57c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {6, 2};
58c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5};
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 8, &points));
619566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes,  triSizes_p2, 2));
629566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(points, triPoints_p2, 8));
63c4762a1bSJed Brown           break;}
64c4762a1bSJed Brown         default:
6563a3b9bcSJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum);
66c4762a1bSJed Brown         }
6730602db0SMatthew G. Knepley       } else if (dim == 2 && simplex && size == 3) {
68c4762a1bSJed Brown         PetscInt triSizes_p3[3]  = {3, 3, 2};
69c4762a1bSJed Brown         PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5};
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(3, &sizes, 8, &points));
729566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(sizes,  triSizes_p3, 3));
739566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(points, triPoints_p3, 8));
7430602db0SMatthew G. Knepley       } else if (dim == 2 && !simplex && size == 2) {
75c4762a1bSJed Brown         PetscInt quadSizes_p2[2]  = {2, 2};
76c4762a1bSJed Brown         PetscInt quadPoints_p2[4] = {2, 3, 0, 1};
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(2, &sizes, 4, &points));
799566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(sizes,  quadSizes_p2, 2));
809566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(points, quadPoints_p2, 4));
81c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
82c4762a1bSJed Brown     }
839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
849566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
859566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
869566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
87c4762a1bSJed Brown   }
889566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(*dm, 0, NULL, &dmDist));
89c4762a1bSJed Brown   if (dmDist) {
909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
91c4762a1bSJed Brown     *dm  = dmDist;
92c4762a1bSJed Brown   }
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh"));
949566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
95c4762a1bSJed Brown   PetscFunctionReturn(0);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user)
99c4762a1bSJed Brown {
100c4762a1bSJed Brown   PetscInt       h, cStart, cEnd, c;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &h));
1049566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd));
105c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
106c4762a1bSJed Brown     /* Could use PetscRand instead */
1079566063dSJacob Faibussowitsch     if (c%2) PetscCall(DMPlexOrientPoint(dm, c, -1));
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown   PetscFunctionReturn(0);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown static PetscErrorCode TestOrientation(DM dm, AppCtx *user)
113c4762a1bSJed Brown {
114c4762a1bSJed Brown   PetscFunctionBeginUser;
1159566063dSJacob Faibussowitsch   PetscCall(ScrambleOrientation(dm, user));
1169566063dSJacob Faibussowitsch   PetscCall(DMPlexOrient(dm));
1179566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-oriented_dm_view"));
118c4762a1bSJed Brown   PetscFunctionReturn(0);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown int main(int argc, char **argv)
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   DM             dm;
124c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
125c4762a1bSJed Brown 
126*327415f7SBarry Smith   PetscFunctionBeginUser;
1279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1289566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1299566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1309566063dSJacob Faibussowitsch   PetscCall(TestOrientation(dm, &user));
1319566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
133b122ec5aSJacob Faibussowitsch   return 0;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*TEST
13730602db0SMatthew G. Knepley   testset:
13830602db0SMatthew G. Knepley     requires: triangle
13930602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view
14030602db0SMatthew G. Knepley 
141c4762a1bSJed Brown     test:
142c4762a1bSJed Brown       suffix: 0
14330602db0SMatthew G. Knepley       args: -test_partition 0
144c4762a1bSJed Brown     test:
145c4762a1bSJed Brown       suffix: 1
146c4762a1bSJed Brown       nsize: 2
147c4762a1bSJed Brown     test:
148c4762a1bSJed Brown       suffix: 2
149c4762a1bSJed Brown       nsize: 2
15030602db0SMatthew G. Knepley       args: -test_num 1
151c4762a1bSJed Brown     test:
152c4762a1bSJed Brown       suffix: 3
153c4762a1bSJed Brown       nsize: 3
15430602db0SMatthew G. Knepley       args: -orientation_view_synchronized
155c4762a1bSJed Brown 
156c4762a1bSJed Brown TEST*/
157