xref: /petsc/src/dm/impls/plex/tests/ex13.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   PetscFunctionBeginUser;
16c4762a1bSJed Brown   options->testPartition = PETSC_TRUE;
17c4762a1bSJed Brown   options->testNum       = 0;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL,0));
221e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
23c4762a1bSJed Brown   PetscFunctionReturn(0);
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
26c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   DM             dmDist = NULL;
2930602db0SMatthew G. Knepley   PetscBool      simplex;
3030602db0SMatthew G. Knepley   PetscInt       dim;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBeginUser;
335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(*dm, &dim));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(*dm, &simplex));
39c4762a1bSJed Brown   if (user->testPartition) {
40c4762a1bSJed Brown     PetscPartitioner part;
41c4762a1bSJed Brown     PetscInt        *sizes  = NULL;
42c4762a1bSJed Brown     PetscInt        *points = NULL;
43c4762a1bSJed Brown     PetscMPIInt      rank, size;
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(comm, &rank));
465f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_size(comm, &size));
47dd400576SPatrick Sanan     if (rank == 0) {
4830602db0SMatthew G. Knepley       if (dim == 2 && simplex && size == 2) {
49c4762a1bSJed Brown         switch (user->testNum) {
50c4762a1bSJed Brown         case 0: {
51c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {4, 4};
52c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4};
53c4762a1bSJed Brown 
545f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscMalloc2(2, &sizes, 8, &points));
555f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscArraycpy(sizes,  triSizes_p2, 2));
565f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscArraycpy(points, triPoints_p2, 8));
57c4762a1bSJed Brown           break;}
58c4762a1bSJed Brown         case 1: {
59c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {6, 2};
60c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5};
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscMalloc2(2, &sizes, 8, &points));
635f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscArraycpy(sizes,  triSizes_p2, 2));
645f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscArraycpy(points, triPoints_p2, 8));
65c4762a1bSJed Brown           break;}
66c4762a1bSJed Brown         default:
6798921bdaSJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum);
68c4762a1bSJed Brown         }
6930602db0SMatthew G. Knepley       } else if (dim == 2 && simplex && size == 3) {
70c4762a1bSJed Brown         PetscInt triSizes_p3[3]  = {3, 3, 2};
71c4762a1bSJed Brown         PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5};
72c4762a1bSJed Brown 
735f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc2(3, &sizes, 8, &points));
745f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(sizes,  triSizes_p3, 3));
755f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(points, triPoints_p3, 8));
7630602db0SMatthew G. Knepley       } else if (dim == 2 && !simplex && size == 2) {
77c4762a1bSJed Brown         PetscInt quadSizes_p2[2]  = {2, 2};
78c4762a1bSJed Brown         PetscInt quadPoints_p2[4] = {2, 3, 0, 1};
79c4762a1bSJed Brown 
805f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc2(2, &sizes, 4, &points));
815f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(sizes,  quadSizes_p2, 2));
825f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(points, quadPoints_p2, 4));
83c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
84c4762a1bSJed Brown     }
855f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(*dm, &part));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerShellSetPartition(part, size, sizes, points));
885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(sizes, points));
89c4762a1bSJed Brown   }
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(*dm, 0, NULL, &dmDist));
91c4762a1bSJed Brown   if (dmDist) {
925f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(dm));
93c4762a1bSJed Brown     *dm  = dmDist;
94c4762a1bSJed Brown   }
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh"));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
97c4762a1bSJed Brown   PetscFunctionReturn(0);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user)
101c4762a1bSJed Brown {
102c4762a1bSJed Brown   PetscInt       h, cStart, cEnd, c;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   PetscFunctionBeginUser;
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetVTKCellHeight(dm, &h));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, h, &cStart, &cEnd));
107c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
108c4762a1bSJed Brown     /* Could use PetscRand instead */
1095f80ce2aSJacob Faibussowitsch     if (c%2) CHKERRQ(DMPlexOrientPoint(dm, c, -1));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown   PetscFunctionReturn(0);
112c4762a1bSJed Brown }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown static PetscErrorCode TestOrientation(DM dm, AppCtx *user)
115c4762a1bSJed Brown {
116c4762a1bSJed Brown   PetscFunctionBeginUser;
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(ScrambleOrientation(dm, user));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexOrient(dm));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-oriented_dm_view"));
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown int main(int argc, char **argv)
124c4762a1bSJed Brown {
125c4762a1bSJed Brown   DM             dm;
126c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
127c4762a1bSJed Brown 
128*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(TestOrientation(dm, &user));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
133*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
134*b122ec5aSJacob Faibussowitsch   return 0;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*TEST
13830602db0SMatthew G. Knepley   testset:
13930602db0SMatthew G. Knepley     requires: triangle
14030602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view
14130602db0SMatthew G. Knepley 
142c4762a1bSJed Brown     test:
143c4762a1bSJed Brown       suffix: 0
14430602db0SMatthew G. Knepley       args: -test_partition 0
145c4762a1bSJed Brown     test:
146c4762a1bSJed Brown       suffix: 1
147c4762a1bSJed Brown       nsize: 2
148c4762a1bSJed Brown     test:
149c4762a1bSJed Brown       suffix: 2
150c4762a1bSJed Brown       nsize: 2
15130602db0SMatthew G. Knepley       args: -test_num 1
152c4762a1bSJed Brown     test:
153c4762a1bSJed Brown       suffix: 3
154c4762a1bSJed Brown       nsize: 3
15530602db0SMatthew G. Knepley       args: -orientation_view_synchronized
156c4762a1bSJed Brown 
157c4762a1bSJed Brown TEST*/
158