xref: /petsc/src/dm/impls/plex/tests/ex11.c (revision 589a23caa660d2a5f330cc8d1ed213e9cfaf51a7)
1c4762a1bSJed Brown static char help[] = "Tests for DMLabel\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static PetscErrorCode TestInsertion()
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   DMLabel        label, label2;
8c4762a1bSJed Brown   const PetscInt values[5] = {0, 3, 4, -1, 176}, N = 10000;
9c4762a1bSJed Brown   PetscInt       i, v;
10c4762a1bSJed Brown   PetscErrorCode ierr;
11c4762a1bSJed Brown 
12c4762a1bSJed Brown   PetscFunctionBegin;
13c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, "Test Label", &label);CHKERRQ(ierr);
14c4762a1bSJed Brown   ierr = DMLabelSetDefaultValue(label, -100);CHKERRQ(ierr);
15c4762a1bSJed Brown   for (i = 0; i < N; ++i) {
16c4762a1bSJed Brown     ierr = DMLabelSetValue(label, i, values[i%5]);CHKERRQ(ierr);
17c4762a1bSJed Brown   }
18c4762a1bSJed Brown   /* Test get in hash mode */
19c4762a1bSJed Brown   for (i = 0; i < N; ++i) {
20c4762a1bSJed Brown     PetscInt val;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     ierr = DMLabelGetValue(label, i, &val);CHKERRQ(ierr);
23c4762a1bSJed Brown     if (val != values[i%5]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %d for point %d should be %d", val, i, values[i%5]);
24c4762a1bSJed Brown   }
25c4762a1bSJed Brown   /* Test stratum */
26c4762a1bSJed Brown   for (v = 0; v < 5; ++v) {
27c4762a1bSJed Brown     IS              stratum;
28c4762a1bSJed Brown     const PetscInt *points;
29c4762a1bSJed Brown     PetscInt        n;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     ierr = DMLabelGetStratumIS(label, values[v], &stratum);CHKERRQ(ierr);
32c4762a1bSJed Brown     if (!stratum) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Stratum %d is empty!", v);
33c4762a1bSJed Brown     ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr);
34c4762a1bSJed Brown     ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr);
35c4762a1bSJed Brown     for (i = 0; i < n; ++i) {
36c4762a1bSJed Brown       if (points[i] != i*5+v) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d should be %d", points[i], i*5+v);
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown     ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr);
39c4762a1bSJed Brown     ierr = ISDestroy(&stratum);CHKERRQ(ierr);
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown   /* Test get in array mode */
42c4762a1bSJed Brown   for (i = 0; i < N; ++i) {
43c4762a1bSJed Brown     PetscInt val;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown     ierr = DMLabelGetValue(label, i, &val);CHKERRQ(ierr);
46c4762a1bSJed Brown     if (val != values[i%5]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %d should be %d", val, values[i%5]);
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   /* Test Duplicate */
49c4762a1bSJed Brown   ierr = DMLabelDuplicate(label, &label2);CHKERRQ(ierr);
50c4762a1bSJed Brown   for (i = 0; i < N; ++i) {
51c4762a1bSJed Brown     PetscInt val;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown     ierr = DMLabelGetValue(label2, i, &val);CHKERRQ(ierr);
54c4762a1bSJed Brown     if (val != values[i%5]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %d should be %d", val, values[i%5]);
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown   ierr = DMLabelDestroy(&label2);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
58c4762a1bSJed Brown   PetscFunctionReturn(0);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown static PetscErrorCode TestEmptyStrata(MPI_Comm comm)
62c4762a1bSJed Brown {
63c4762a1bSJed Brown   DM               dm, dmDist;
64c4762a1bSJed Brown   PetscPartitioner part;
65c4762a1bSJed Brown   PetscInt         c0[6]  = {2,3,6,7,9,11};
66c4762a1bSJed Brown   PetscInt         c1[6]  = {4,5,7,8,10,12};
67c4762a1bSJed Brown   PetscInt         c2[4]  = {13,15,19,21};
68c4762a1bSJed Brown   PetscInt         c3[4]  = {14,16,20,22};
69c4762a1bSJed Brown   PetscInt         c4[4]  = {15,17,21,23};
70c4762a1bSJed Brown   PetscInt         c5[4]  = {16,18,22,24};
71c4762a1bSJed Brown   PetscInt         c6[4]  = {13,14,19,20};
72c4762a1bSJed Brown   PetscInt         c7[4]  = {15,16,21,22};
73c4762a1bSJed Brown   PetscInt         c8[4]  = {17,18,23,24};
74c4762a1bSJed Brown   PetscInt         c9[4]  = {13,14,15,16};
75c4762a1bSJed Brown   PetscInt         c10[4] = {15,16,17,18};
76c4762a1bSJed Brown   PetscInt         c11[4] = {19,20,21,22};
77c4762a1bSJed Brown   PetscInt         c12[4] = {21,22,23,24};
78c4762a1bSJed Brown   PetscInt         dim    = 3;
79c4762a1bSJed Brown   PetscMPIInt      rank;
80c4762a1bSJed Brown   PetscErrorCode   ierr;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBegin;
83c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
84c4762a1bSJed Brown   /* A 3D box with two adjacent cells, sharing one face and four vertices */
85c4762a1bSJed Brown   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
86c4762a1bSJed Brown   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
88c4762a1bSJed Brown   if (!rank) {
89c4762a1bSJed Brown     ierr = DMPlexSetChart(dm, 0, 25);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 0, 6);CHKERRQ(ierr);
91c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 1, 6);CHKERRQ(ierr);
92c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 2, 4);CHKERRQ(ierr);
93c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 3, 4);CHKERRQ(ierr);
94c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 4, 4);CHKERRQ(ierr);
95c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 5, 4);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 6, 4);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 7, 4);CHKERRQ(ierr);
98c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 8, 4);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 9, 4);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 10, 4);CHKERRQ(ierr);
101c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 11, 4);CHKERRQ(ierr);
102c4762a1bSJed Brown     ierr = DMPlexSetConeSize(dm, 12, 4);CHKERRQ(ierr);
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
105c4762a1bSJed Brown   if (!rank) {
106c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 0, c0);CHKERRQ(ierr);
107c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 1, c1);CHKERRQ(ierr);
108c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 2, c2);CHKERRQ(ierr);
109c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 3, c3);CHKERRQ(ierr);
110c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 4, c4);CHKERRQ(ierr);
111c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 5, c5);CHKERRQ(ierr);
112c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 6, c6);CHKERRQ(ierr);
113c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 7, c7);CHKERRQ(ierr);
114c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 8, c8);CHKERRQ(ierr);
115c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 9, c9);CHKERRQ(ierr);
116c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 10, c10);CHKERRQ(ierr);
117c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 11, c11);CHKERRQ(ierr);
118c4762a1bSJed Brown     ierr = DMPlexSetCone(dm, 12, c12);CHKERRQ(ierr);
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown   ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);
121c4762a1bSJed Brown   /* Create a user managed depth label, so that we can leave out edges */
122c4762a1bSJed Brown   {
123c4762a1bSJed Brown     DMLabel label;
124c4762a1bSJed Brown     PetscInt numValues, maxValues = 0, v;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     ierr = DMCreateLabel(dm, "depth");CHKERRQ(ierr);
127c4762a1bSJed Brown     ierr = DMPlexGetDepthLabel(dm, &label);CHKERRQ(ierr);
128c4762a1bSJed Brown     if (!rank) {
129c4762a1bSJed Brown       PetscInt i;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown       for (i = 0; i < 25; ++i) {
132c4762a1bSJed Brown         if (i < 2)       {ierr = DMLabelSetValue(label, i, 3);CHKERRQ(ierr);}
133c4762a1bSJed Brown         else if (i < 13) {ierr = DMLabelSetValue(label, i, 2);CHKERRQ(ierr);}
134c4762a1bSJed Brown         else             {
135c4762a1bSJed Brown           if (i==13) {ierr = DMLabelAddStratum(label, 1);CHKERRQ(ierr);}
136c4762a1bSJed Brown           ierr = DMLabelSetValue(label, i, 0);CHKERRQ(ierr);
137c4762a1bSJed Brown         }
138c4762a1bSJed Brown       }
139c4762a1bSJed Brown     }
140c4762a1bSJed Brown     ierr = DMLabelGetNumValues(label, &numValues);CHKERRQ(ierr);
141c4762a1bSJed Brown     ierr = MPI_Allreduce(&numValues, &maxValues, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
142c4762a1bSJed Brown     for (v = numValues; v < maxValues; ++v) {ierr = DMLabelAddStratum(label,v);CHKERRQ(ierr);}
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown   {
145c4762a1bSJed Brown     DMLabel label;
146c4762a1bSJed Brown     ierr = DMPlexGetDepthLabel(dm, &label);CHKERRQ(ierr);
147c4762a1bSJed Brown     ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown   ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = DMPlexDistribute(dm, 1, NULL, &dmDist);CHKERRQ(ierr);
152c4762a1bSJed Brown   if (dmDist) {
153c4762a1bSJed Brown     ierr = DMDestroy(&dm);CHKERRQ(ierr);
154c4762a1bSJed Brown     dm   = dmDist;
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown   {
157c4762a1bSJed Brown     DMLabel label;
158c4762a1bSJed Brown     ierr = DMPlexGetDepthLabel(dm, &label);CHKERRQ(ierr);
159c4762a1bSJed Brown     ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown   /* Create a cell vector */
162c4762a1bSJed Brown   {
163c4762a1bSJed Brown     Vec          v;
164c4762a1bSJed Brown     PetscSection s;
165c4762a1bSJed Brown     PetscInt     numComp[] = {1};
166c4762a1bSJed Brown     PetscInt     dof[]     = {0,0,0,1};
167c4762a1bSJed Brown     PetscInt     N;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown     ierr = DMSetNumFields(dm, 1);CHKERRQ(ierr);
170c4762a1bSJed Brown     ierr = DMPlexCreateSection(dm, NULL, numComp, dof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr);
171c4762a1bSJed Brown     ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr);
172c4762a1bSJed Brown     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
173c4762a1bSJed Brown     ierr = DMCreateGlobalVector(dm, &v);CHKERRQ(ierr);
174c4762a1bSJed Brown     ierr = VecGetSize(v, &N);CHKERRQ(ierr);
175c4762a1bSJed Brown     if (N != 2) {
176c4762a1bSJed Brown       ierr = DMView(dm, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
177c4762a1bSJed Brown       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "FAIL: Vector size %d != 2\n", N);
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown static PetscErrorCode TestDistribution(MPI_Comm comm)
186c4762a1bSJed Brown {
187c4762a1bSJed Brown   DM               dm, dmDist;
188c4762a1bSJed Brown   PetscPartitioner part;
189c4762a1bSJed Brown   DMLabel          label;
190c4762a1bSJed Brown   char             filename[2048];
191c4762a1bSJed Brown   const char      *name    = "test label";
192c4762a1bSJed Brown   PetscInt         overlap = 0, cStart, cEnd, c;
193c4762a1bSJed Brown   PetscMPIInt      rank;
194c4762a1bSJed Brown   PetscBool        flg;
195c4762a1bSJed Brown   PetscErrorCode   ierr;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBegin;
198c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
199*589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL, NULL, "-filename", filename, sizeof(filename), &flg);CHKERRQ(ierr);
200c4762a1bSJed Brown   if (!flg) PetscFunctionReturn(0);
201c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL, NULL, "-overlap", &overlap, NULL);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm); CHKERRQ(ierr);
203c4762a1bSJed Brown   ierr = DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE);CHKERRQ(ierr);
204c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
205c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
207c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
208c4762a1bSJed Brown     ierr = DMLabelSetValue(label, c, c);CHKERRQ(ierr);
209c4762a1bSJed Brown   }
210c4762a1bSJed Brown   ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
211c4762a1bSJed Brown   ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = DMPlexDistribute(dm, overlap, NULL, &dmDist);CHKERRQ(ierr);
214c4762a1bSJed Brown   if (dmDist) {
215c4762a1bSJed Brown     ierr = DMDestroy(&dm);CHKERRQ(ierr);
216c4762a1bSJed Brown     dm   = dmDist;
217c4762a1bSJed Brown   }
218c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
220c4762a1bSJed Brown   ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
221c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
222c4762a1bSJed Brown   PetscFunctionReturn(0);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown int main(int argc, char **argv)
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   PetscErrorCode ierr;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
230c4762a1bSJed Brown   /*ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);*/
231c4762a1bSJed Brown   ierr = TestInsertion();CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = TestEmptyStrata(PETSC_COMM_WORLD);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = TestDistribution(PETSC_COMM_WORLD);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = PetscFinalize();
235c4762a1bSJed Brown   return ierr;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown /*TEST
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   test:
241c4762a1bSJed Brown     suffix: 0
242c4762a1bSJed Brown   test:
243c4762a1bSJed Brown     suffix: 1
244c4762a1bSJed Brown     nsize: 2
245c4762a1bSJed Brown     args: -petscpartitioner_type simple
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   testset:
248c4762a1bSJed Brown     suffix: gmsh
249c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -petscpartitioner_type simple
250c4762a1bSJed Brown     test:
251c4762a1bSJed Brown       suffix: 1
252c4762a1bSJed Brown       nsize: 1
253c4762a1bSJed Brown     test:
254c4762a1bSJed Brown       suffix: 2
255c4762a1bSJed Brown       nsize: 2
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   testset:
258c4762a1bSJed Brown     suffix: exodusii
259c4762a1bSJed Brown     requires: exodusii
260c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/2Dgrd.exo -petscpartitioner_type simple
261c4762a1bSJed Brown     test:
262c4762a1bSJed Brown       suffix: 1
263c4762a1bSJed Brown       nsize: 1
264c4762a1bSJed Brown     test:
265c4762a1bSJed Brown       suffix: 2
266c4762a1bSJed Brown       nsize: 2
267c4762a1bSJed Brown 
268c4762a1bSJed Brown TEST*/
269