xref: /petsc/src/dm/impls/plex/tests/ex13.c (revision 4d81c09ef4cefe32a382d98cf0a16804df0ea70d)
1 static char help[] = "Orient a mesh in parallel\n\n";
2 
3 #include <petscdmplex.h>
4 
5 typedef struct {
6   /* Domain and mesh definition */
7   PetscBool testPartition; /* Use a fixed partitioning for testing */
8   PetscInt  testNum;       /* Labels the different test partitions */
9 } AppCtx;
10 
11 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12 {
13   PetscErrorCode ierr;
14 
15   PetscFunctionBeginUser;
16   options->testPartition = PETSC_TRUE;
17   options->testNum       = 0;
18 
19   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
20   ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr);
21   ierr = PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr);
22   ierr = PetscOptionsEnd();CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
27 {
28   DM             dmDist = NULL;
29   PetscBool      simplex;
30   PetscInt       dim;
31   PetscErrorCode ierr;
32 
33   PetscFunctionBeginUser;
34   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
35   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
36   ierr = DMPlexDistributeSetDefault(*dm, PETSC_FALSE);CHKERRQ(ierr);
37   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
38   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
39   ierr = DMPlexIsSimplex(*dm, &simplex);CHKERRQ(ierr);
40   if (user->testPartition) {
41     PetscPartitioner part;
42     PetscInt        *sizes  = NULL;
43     PetscInt        *points = NULL;
44     PetscMPIInt      rank, size;
45 
46     ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
47     ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
48     if (rank == 0) {
49       if (dim == 2 && simplex && size == 2) {
50         switch (user->testNum) {
51         case 0: {
52           PetscInt triSizes_p2[2]  = {4, 4};
53           PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4};
54 
55           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
56           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
57           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
58           break;}
59         case 1: {
60           PetscInt triSizes_p2[2]  = {6, 2};
61           PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5};
62 
63           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
64           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
65           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
66           break;}
67         default:
68           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum);
69         }
70       } else if (dim == 2 && simplex && size == 3) {
71         PetscInt triSizes_p3[3]  = {3, 3, 2};
72         PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5};
73 
74         ierr = PetscMalloc2(3, &sizes, 8, &points);CHKERRQ(ierr);
75         ierr = PetscArraycpy(sizes,  triSizes_p3, 3);CHKERRQ(ierr);
76         ierr = PetscArraycpy(points, triPoints_p3, 8);CHKERRQ(ierr);
77       } else if (dim == 2 && !simplex && size == 2) {
78         PetscInt quadSizes_p2[2]  = {2, 2};
79         PetscInt quadPoints_p2[4] = {2, 3, 0, 1};
80 
81         ierr = PetscMalloc2(2, &sizes, 4, &points);CHKERRQ(ierr);
82         ierr = PetscArraycpy(sizes,  quadSizes_p2, 2);CHKERRQ(ierr);
83         ierr = PetscArraycpy(points, quadPoints_p2, 4);CHKERRQ(ierr);
84       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
85     }
86     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
87     ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
88     ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
89     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
90   }
91   ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
92   if (dmDist) {
93     ierr = DMDestroy(dm);CHKERRQ(ierr);
94     *dm  = dmDist;
95   }
96   ierr = PetscObjectSetName((PetscObject) *dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh");CHKERRQ(ierr);
97   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user)
102 {
103   PetscInt       h, cStart, cEnd, c;
104   PetscErrorCode ierr;
105 
106   PetscFunctionBeginUser;
107   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
108   ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr);
109   for (c = cStart; c < cEnd; ++c) {
110     /* Could use PetscRand instead */
111     if (c%2) {ierr = DMPlexOrientPoint(dm, c, -1);CHKERRQ(ierr);}
112   }
113   PetscFunctionReturn(0);
114 }
115 
116 static PetscErrorCode TestOrientation(DM dm, AppCtx *user)
117 {
118   PetscErrorCode ierr;
119 
120   PetscFunctionBeginUser;
121   ierr = ScrambleOrientation(dm, user);CHKERRQ(ierr);
122   ierr = DMPlexOrient(dm);CHKERRQ(ierr);
123   ierr = DMViewFromOptions(dm, NULL, "-oriented_dm_view");CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 int main(int argc, char **argv)
128 {
129   DM             dm;
130   AppCtx         user; /* user-defined work context */
131   PetscErrorCode ierr;
132 
133   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
134   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
135   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
136   ierr = TestOrientation(dm, &user);CHKERRQ(ierr);
137   ierr = DMDestroy(&dm);CHKERRQ(ierr);
138   ierr = PetscFinalize();
139   return ierr;
140 }
141 
142 /*TEST
143   testset:
144     requires: triangle
145     args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view
146 
147     test:
148       suffix: 0
149       args: -test_partition 0
150     test:
151       suffix: 1
152       nsize: 2
153     test:
154       suffix: 2
155       nsize: 2
156       args: -test_num 1
157     test:
158       suffix: 3
159       nsize: 3
160       args: -orientation_view_synchronized
161 
162 TEST*/
163