xref: /petsc/src/dm/impls/plex/tests/ex13.c (revision 2fa40bb9206b96114faa7cb222621ec184d31cd2)
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 = DMSetFromOptions(*dm);CHKERRQ(ierr);
37   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
38   ierr = DMPlexIsSimplex(*dm, &simplex);CHKERRQ(ierr);
39   if (user->testPartition) {
40     PetscPartitioner part;
41     PetscInt        *sizes  = NULL;
42     PetscInt        *points = NULL;
43     PetscMPIInt      rank, size;
44 
45     ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
46     ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
47     if (rank == 0) {
48       if (dim == 2 && simplex && size == 2) {
49         switch (user->testNum) {
50         case 0: {
51           PetscInt triSizes_p2[2]  = {4, 4};
52           PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4};
53 
54           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
55           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
56           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
57           break;}
58         case 1: {
59           PetscInt triSizes_p2[2]  = {6, 2};
60           PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5};
61 
62           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
63           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
64           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
65           break;}
66         default:
67           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum);
68         }
69       } else if (dim == 2 && simplex && size == 3) {
70         PetscInt triSizes_p3[3]  = {3, 3, 2};
71         PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5};
72 
73         ierr = PetscMalloc2(3, &sizes, 8, &points);CHKERRQ(ierr);
74         ierr = PetscArraycpy(sizes,  triSizes_p3, 3);CHKERRQ(ierr);
75         ierr = PetscArraycpy(points, triPoints_p3, 8);CHKERRQ(ierr);
76       } else if (dim == 2 && !simplex && size == 2) {
77         PetscInt quadSizes_p2[2]  = {2, 2};
78         PetscInt quadPoints_p2[4] = {2, 3, 0, 1};
79 
80         ierr = PetscMalloc2(2, &sizes, 4, &points);CHKERRQ(ierr);
81         ierr = PetscArraycpy(sizes,  quadSizes_p2, 2);CHKERRQ(ierr);
82         ierr = PetscArraycpy(points, quadPoints_p2, 4);CHKERRQ(ierr);
83       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
84     }
85     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
86     ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
87     ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
88     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
89   }
90   ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
91   if (dmDist) {
92     ierr = DMDestroy(dm);CHKERRQ(ierr);
93     *dm  = dmDist;
94   }
95   ierr = PetscObjectSetName((PetscObject) *dm, simplex ? "Simplicial Mesh" : "Tensor Product Mesh");CHKERRQ(ierr);
96   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 static PetscErrorCode ScrambleOrientation(DM dm, AppCtx *user)
101 {
102   PetscInt       h, cStart, cEnd, c;
103   PetscErrorCode ierr;
104 
105   PetscFunctionBeginUser;
106   ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr);
107   ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr);
108   for (c = cStart; c < cEnd; ++c) {
109     /* Could use PetscRand instead */
110     if (c%2) {ierr = DMPlexOrientPoint(dm, c, -1);CHKERRQ(ierr);}
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 static PetscErrorCode TestOrientation(DM dm, AppCtx *user)
116 {
117   PetscErrorCode ierr;
118 
119   PetscFunctionBeginUser;
120   ierr = ScrambleOrientation(dm, user);CHKERRQ(ierr);
121   ierr = DMPlexOrient(dm);CHKERRQ(ierr);
122   ierr = DMViewFromOptions(dm, NULL, "-oriented_dm_view");CHKERRQ(ierr);
123   PetscFunctionReturn(0);
124 }
125 
126 int main(int argc, char **argv)
127 {
128   DM             dm;
129   AppCtx         user; /* user-defined work context */
130   PetscErrorCode ierr;
131 
132   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
133   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
134   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
135   ierr = TestOrientation(dm, &user);CHKERRQ(ierr);
136   ierr = DMDestroy(&dm);CHKERRQ(ierr);
137   ierr = PetscFinalize();
138   return ierr;
139 }
140 
141 /*TEST
142   testset:
143     requires: triangle
144     args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -oriented_dm_view ascii::ascii_info_detail -orientation_view
145 
146     test:
147       suffix: 0
148       args: -test_partition 0
149     test:
150       suffix: 1
151       nsize: 2
152     test:
153       suffix: 2
154       nsize: 2
155       args: -test_num 1
156     test:
157       suffix: 3
158       nsize: 3
159       args: -orientation_view_synchronized
160 
161 TEST*/
162