xref: /petsc/src/sys/tests/ex62.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c0888a1bSJDBetteridge 
2c0888a1bSJDBetteridge static char help[] = "Tests `PetscGarbageKeySortedIntersect()`\n\n";
3c0888a1bSJDBetteridge 
4c0888a1bSJDBetteridge #include <petscsys.h>
5c0888a1bSJDBetteridge #include <petsc/private/garbagecollector.h>
6c0888a1bSJDBetteridge 
7c0888a1bSJDBetteridge /* This program tests `PetscGarbageKeySortedIntersect(), which is the
8c0888a1bSJDBetteridge    public (MPI) interface to
9c0888a1bSJDBetteridge    `PetscErrorCode GarbageKeySortedIntersect_Private()`.
10c0888a1bSJDBetteridge    Sets are sent packed in arrays, with the first entry as the number of
11c0888a1bSJDBetteridge    set elements and the sets the remaining elements. This is because the
12c0888a1bSJDBetteridge    MPI reduction operation must have the call signature:
13c0888a1bSJDBetteridge    void PetscGarbageKeySortedIntersect(void *inset, void *inoutset, PetscMPIInt *length, MPI_Datatype *dtype)
14c0888a1bSJDBetteridge    This is a thin wrapper for the private routine:
15c0888a1bSJDBetteridge    PetscErrorCode GarbageKeySortedIntersect_Private(PetscInt64 seta[], PetscInt *lena, PetscInt64 setb[], PetscInt lenb)
16c0888a1bSJDBetteridge    Where
17c0888a1bSJDBetteridge    seta = (PetscInt64 *)inoutset;
18c0888a1bSJDBetteridge    setb = (PetscInt64 *)inset;
19c0888a1bSJDBetteridge    And the arguments are passed as:
20c0888a1bSJDBetteridge    &seta[1], (PetscInt *)&seta[0], &setb[1], (PetscInt)setb[0]
21c0888a1bSJDBetteridge */
22c0888a1bSJDBetteridge 
23c0888a1bSJDBetteridge /* Populate a set with upto the first 49 unique Fibonnaci numbers */
24c0888a1bSJDBetteridge PetscErrorCode Fibonnaci(PetscInt64 **set, PetscInt n)
25c0888a1bSJDBetteridge {
26c0888a1bSJDBetteridge   PetscInt   ii;
27c0888a1bSJDBetteridge   PetscInt64 fib[] = {1,        2,        3,        5,        8,         13,        21,        34,        55,        89,         144,        233,        377,        610,        987,        1597,    2584,
28c0888a1bSJDBetteridge                       4181,     6765,     10946,    17711,    28657,     46368,     75025,     121393,    196418,    317811,     514229,     832040,     1346269,    2178309,    3524578,    5702887, 9227465,
29c0888a1bSJDBetteridge                       14930352, 24157817, 39088169, 63245986, 102334155, 165580141, 267914296, 433494437, 701408733, 1134903170, 1836311903, 2971215073, 4807526976, 7778742049, 12586269025};
30c0888a1bSJDBetteridge 
31c0888a1bSJDBetteridge   PetscFunctionBeginUser;
32c0888a1bSJDBetteridge   PetscAssert((n < 50), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "n must be less than 50\n");
33c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
34c0888a1bSJDBetteridge   (*set)[0] = (PetscInt64)n;
35c0888a1bSJDBetteridge   for (ii = 0; ii < n; ii++) { (*set)[ii + 1] = fib[ii]; }
36*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37c0888a1bSJDBetteridge }
38c0888a1bSJDBetteridge 
39c0888a1bSJDBetteridge /* Populate a set with Square numbers */
40c0888a1bSJDBetteridge PetscErrorCode Square(PetscInt64 **set, PetscInt n)
41c0888a1bSJDBetteridge {
42c0888a1bSJDBetteridge   PetscInt64 ii;
43c0888a1bSJDBetteridge 
44c0888a1bSJDBetteridge   PetscFunctionBeginUser;
45c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
46c0888a1bSJDBetteridge   (*set)[0] = (PetscInt64)n;
47c0888a1bSJDBetteridge   for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii; }
48*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49c0888a1bSJDBetteridge }
50c0888a1bSJDBetteridge 
51c0888a1bSJDBetteridge /* Populate a set with Cube numbers */
52c0888a1bSJDBetteridge PetscErrorCode Cube(PetscInt64 **set, PetscInt n)
53c0888a1bSJDBetteridge {
54c0888a1bSJDBetteridge   PetscInt64 ii;
55c0888a1bSJDBetteridge 
56c0888a1bSJDBetteridge   PetscFunctionBeginUser;
57c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
58c0888a1bSJDBetteridge   (*set)[0] = (PetscInt64)n;
59c0888a1bSJDBetteridge   for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii * ii; }
60*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61c0888a1bSJDBetteridge }
62c0888a1bSJDBetteridge 
63c0888a1bSJDBetteridge /* Populate a set with numbers to sixth power */
64c0888a1bSJDBetteridge PetscErrorCode Sixth(PetscInt64 **set, PetscInt n)
65c0888a1bSJDBetteridge {
66c0888a1bSJDBetteridge   PetscInt64 ii;
67c0888a1bSJDBetteridge 
68c0888a1bSJDBetteridge   PetscFunctionBeginUser;
69c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
70c0888a1bSJDBetteridge   (*set)[0] = (PetscInt64)n;
71c0888a1bSJDBetteridge   for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii * ii * ii * ii * ii; }
72*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73c0888a1bSJDBetteridge }
74c0888a1bSJDBetteridge 
75c0888a1bSJDBetteridge /* Print out the contents of a set */
76c0888a1bSJDBetteridge PetscErrorCode PrintSet(PetscInt64 *set)
77c0888a1bSJDBetteridge {
78c0888a1bSJDBetteridge   char     text[64];
79c0888a1bSJDBetteridge   PetscInt ii;
80c0888a1bSJDBetteridge 
81c0888a1bSJDBetteridge   PetscFunctionBeginUser;
82c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "["));
83c0888a1bSJDBetteridge   for (ii = 1; ii <= (PetscInt)set[0]; ii++) {
84c0888a1bSJDBetteridge     PetscCall(PetscFormatConvert(" %" PetscInt64_FMT ",", text));
85c0888a1bSJDBetteridge     PetscCall(PetscPrintf(PETSC_COMM_WORLD, text, set[ii]));
86c0888a1bSJDBetteridge   }
87c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
88*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89c0888a1bSJDBetteridge }
90c0888a1bSJDBetteridge 
91c0888a1bSJDBetteridge /* Check set equality */
92c0888a1bSJDBetteridge PetscErrorCode AssertSetsEqual(PetscInt64 *set, PetscInt64 *true_set)
93c0888a1bSJDBetteridge {
94c0888a1bSJDBetteridge   PetscInt ii;
95c0888a1bSJDBetteridge 
96c0888a1bSJDBetteridge   PetscFunctionBeginUser;
97c0888a1bSJDBetteridge   PetscAssert((set[0] == true_set[0]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets of different sizes");
98c0888a1bSJDBetteridge   for (ii = 1; ii < set[0] + 1; ii++) { PetscAssert((set[ii] == true_set[ii]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets are different"); }
99*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100c0888a1bSJDBetteridge }
101c0888a1bSJDBetteridge 
102c0888a1bSJDBetteridge /* Tests functionality when two enpty sets are passed */
103c0888a1bSJDBetteridge PetscErrorCode test_empty_empty()
104c0888a1bSJDBetteridge {
105c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
106c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
107c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
108c0888a1bSJDBetteridge 
109c0888a1bSJDBetteridge   PetscFunctionBeginUser;
110c0888a1bSJDBetteridge 
111c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
112c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
113c0888a1bSJDBetteridge 
114c0888a1bSJDBetteridge   set_a[0] = 0;
115c0888a1bSJDBetteridge 
116c0888a1bSJDBetteridge   set_b[0] = 0;
117c0888a1bSJDBetteridge 
118c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
119c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
120c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
121c0888a1bSJDBetteridge 
122c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
123c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
124c0888a1bSJDBetteridge 
125*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126c0888a1bSJDBetteridge }
127c0888a1bSJDBetteridge 
128c0888a1bSJDBetteridge /* Tests functionality when seta is empty */
129c0888a1bSJDBetteridge PetscErrorCode test_a_empty()
130c0888a1bSJDBetteridge {
131c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
132c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
133c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
134c0888a1bSJDBetteridge 
135c0888a1bSJDBetteridge   PetscFunctionBeginUser;
136c0888a1bSJDBetteridge 
137c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
138c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_b));
139c0888a1bSJDBetteridge 
140c0888a1bSJDBetteridge   set_a[0] = 0;
141c0888a1bSJDBetteridge 
142c0888a1bSJDBetteridge   set_b[0] = 1;
143c0888a1bSJDBetteridge   set_b[1] = 1;
144c0888a1bSJDBetteridge 
145c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
146c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
147c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
148c0888a1bSJDBetteridge 
149c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
150c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
151c0888a1bSJDBetteridge 
152*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153c0888a1bSJDBetteridge }
154c0888a1bSJDBetteridge 
155c0888a1bSJDBetteridge /* Tests functionality when setb is empty */
156c0888a1bSJDBetteridge PetscErrorCode test_b_empty()
157c0888a1bSJDBetteridge {
158c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
159c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
160c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
161c0888a1bSJDBetteridge 
162c0888a1bSJDBetteridge   PetscFunctionBeginUser;
163c0888a1bSJDBetteridge 
164c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_a));
165c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
166c0888a1bSJDBetteridge 
167c0888a1bSJDBetteridge   set_a[0] = 1;
168c0888a1bSJDBetteridge   set_a[1] = 1;
169c0888a1bSJDBetteridge 
170c0888a1bSJDBetteridge   set_b[0] = 0;
171c0888a1bSJDBetteridge 
172c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
173c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
174c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
175c0888a1bSJDBetteridge 
176c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
177c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
178c0888a1bSJDBetteridge 
179*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180c0888a1bSJDBetteridge }
181c0888a1bSJDBetteridge 
182c0888a1bSJDBetteridge /* Tests functionality when both sets are identical */
183c0888a1bSJDBetteridge PetscErrorCode test_identical()
184c0888a1bSJDBetteridge {
185c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
186c0888a1bSJDBetteridge   PetscInt64  truth[] = {3, 1, 4, 9};
187c0888a1bSJDBetteridge   PetscMPIInt length  = 4;
188c0888a1bSJDBetteridge 
189c0888a1bSJDBetteridge   PetscFunctionBeginUser;
190c0888a1bSJDBetteridge 
191c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
192c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_b));
193c0888a1bSJDBetteridge 
194c0888a1bSJDBetteridge   set_a[0] = 3;
195c0888a1bSJDBetteridge   set_a[1] = 1;
196c0888a1bSJDBetteridge   set_a[2] = 4;
197c0888a1bSJDBetteridge   set_a[3] = 9;
198c0888a1bSJDBetteridge 
199c0888a1bSJDBetteridge   set_b[0] = 3;
200c0888a1bSJDBetteridge   set_b[1] = 1;
201c0888a1bSJDBetteridge   set_b[2] = 4;
202c0888a1bSJDBetteridge   set_b[3] = 9;
203c0888a1bSJDBetteridge 
204c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
205c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
206c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
207c0888a1bSJDBetteridge 
208c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
209c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
210c0888a1bSJDBetteridge 
211*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212c0888a1bSJDBetteridge }
213c0888a1bSJDBetteridge 
214c0888a1bSJDBetteridge /* Tests functionality when sets have no elements in common */
215c0888a1bSJDBetteridge PetscErrorCode test_disjoint()
216c0888a1bSJDBetteridge {
217c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
218c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
219c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
220c0888a1bSJDBetteridge 
221c0888a1bSJDBetteridge   PetscFunctionBeginUser;
222c0888a1bSJDBetteridge 
223c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
224c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_b));
225c0888a1bSJDBetteridge 
226c0888a1bSJDBetteridge   set_a[0] = 3;
227c0888a1bSJDBetteridge   set_a[1] = 1;
228c0888a1bSJDBetteridge   set_a[2] = 4;
229c0888a1bSJDBetteridge   set_a[3] = 9;
230c0888a1bSJDBetteridge 
231c0888a1bSJDBetteridge   set_b[0] = 3;
232c0888a1bSJDBetteridge   set_b[1] = 2;
233c0888a1bSJDBetteridge   set_b[2] = 6;
234c0888a1bSJDBetteridge   set_b[3] = 8;
235c0888a1bSJDBetteridge 
236c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
237c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
238c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
239c0888a1bSJDBetteridge 
240c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
241c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
242c0888a1bSJDBetteridge 
243*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244c0888a1bSJDBetteridge }
245c0888a1bSJDBetteridge 
246c0888a1bSJDBetteridge /* Tests functionality when sets only have one element in common */
247c0888a1bSJDBetteridge PetscErrorCode test_single_common()
248c0888a1bSJDBetteridge {
249c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
250c0888a1bSJDBetteridge   PetscInt64  truth[] = {1, 4};
251c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
252c0888a1bSJDBetteridge 
253c0888a1bSJDBetteridge   PetscFunctionBeginUser;
254c0888a1bSJDBetteridge 
255c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
256c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_b));
257c0888a1bSJDBetteridge 
258c0888a1bSJDBetteridge   set_a[0] = 3;
259c0888a1bSJDBetteridge   set_a[1] = 1;
260c0888a1bSJDBetteridge   set_a[2] = 4;
261c0888a1bSJDBetteridge   set_a[3] = 9;
262c0888a1bSJDBetteridge 
263c0888a1bSJDBetteridge   set_b[0] = 3;
264c0888a1bSJDBetteridge   set_b[1] = 2;
265c0888a1bSJDBetteridge   set_b[2] = 4;
266c0888a1bSJDBetteridge   set_b[3] = 6;
267c0888a1bSJDBetteridge   set_b[4] = 8;
268c0888a1bSJDBetteridge 
269c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
270c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
271c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
272c0888a1bSJDBetteridge 
273c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
274c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
275c0888a1bSJDBetteridge 
276*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277c0888a1bSJDBetteridge }
278c0888a1bSJDBetteridge 
279c0888a1bSJDBetteridge /* Specific test case flagged by PETSc issue #1247 */
280c0888a1bSJDBetteridge PetscErrorCode test_issue_1247()
281c0888a1bSJDBetteridge {
282c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
283c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
284c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
285c0888a1bSJDBetteridge 
286c0888a1bSJDBetteridge   PetscFunctionBeginUser;
287c0888a1bSJDBetteridge 
288c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(3, &set_a));
289c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_b));
290c0888a1bSJDBetteridge 
291c0888a1bSJDBetteridge   set_a[0] = 2;
292c0888a1bSJDBetteridge   set_a[1] = 2;
293c0888a1bSJDBetteridge   set_a[2] = 3;
294c0888a1bSJDBetteridge 
295c0888a1bSJDBetteridge   set_b[0] = 1;
296c0888a1bSJDBetteridge   set_b[1] = 1;
297c0888a1bSJDBetteridge 
298c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
299c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
300c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
301c0888a1bSJDBetteridge 
302c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
303c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
304c0888a1bSJDBetteridge 
305*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306c0888a1bSJDBetteridge }
307c0888a1bSJDBetteridge 
308c0888a1bSJDBetteridge /* Tests functionality when seta is empty and setb is large */
309c0888a1bSJDBetteridge PetscErrorCode test_empty_big()
310c0888a1bSJDBetteridge {
311c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
312c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
313c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
314c0888a1bSJDBetteridge 
315c0888a1bSJDBetteridge   PetscFunctionBeginUser;
316c0888a1bSJDBetteridge 
317c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
318c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
319c0888a1bSJDBetteridge 
320c0888a1bSJDBetteridge   set_a[0] = 0;
321c0888a1bSJDBetteridge 
322c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
323c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
324c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
325c0888a1bSJDBetteridge 
326c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
327c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
328c0888a1bSJDBetteridge 
329*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
330c0888a1bSJDBetteridge }
331c0888a1bSJDBetteridge 
332c0888a1bSJDBetteridge /* Tests functionality when seta is small and setb is large */
333c0888a1bSJDBetteridge PetscErrorCode test_small_big()
334c0888a1bSJDBetteridge {
335c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
336c0888a1bSJDBetteridge   PetscInt64  truth[] = {3, 1, 4, 9};
337c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
338c0888a1bSJDBetteridge 
339c0888a1bSJDBetteridge   PetscFunctionBeginUser;
340c0888a1bSJDBetteridge 
341c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_a));
342c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
343c0888a1bSJDBetteridge 
344c0888a1bSJDBetteridge   set_a[0] = 4;
345c0888a1bSJDBetteridge   set_a[1] = 1;
346c0888a1bSJDBetteridge   set_a[2] = 4;
347c0888a1bSJDBetteridge   set_a[3] = 8;
348c0888a1bSJDBetteridge   set_a[4] = 9;
349c0888a1bSJDBetteridge 
350c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
351c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
352c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
353c0888a1bSJDBetteridge 
354c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
355c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
356c0888a1bSJDBetteridge 
357*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358c0888a1bSJDBetteridge }
359c0888a1bSJDBetteridge 
360c0888a1bSJDBetteridge /* Tests functionality when seta is medium sized and setb is large */
361c0888a1bSJDBetteridge PetscErrorCode test_moderate_big()
362c0888a1bSJDBetteridge {
363c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
364c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 144};
365c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
366c0888a1bSJDBetteridge 
367c0888a1bSJDBetteridge   PetscFunctionBeginUser;
368c0888a1bSJDBetteridge 
369c0888a1bSJDBetteridge   PetscCall(Fibonnaci(&set_a, 49));
370c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
371c0888a1bSJDBetteridge 
372c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
373c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
374c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
375c0888a1bSJDBetteridge 
376c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
377c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
378c0888a1bSJDBetteridge 
379*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
380c0888a1bSJDBetteridge }
381c0888a1bSJDBetteridge 
382c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large */
383c0888a1bSJDBetteridge PetscErrorCode test_big_big()
384c0888a1bSJDBetteridge {
385c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
386c0888a1bSJDBetteridge   PetscInt64 *truth;
387c0888a1bSJDBetteridge   PetscMPIInt length = 1;
388c0888a1bSJDBetteridge 
389c0888a1bSJDBetteridge   PetscFunctionBeginUser;
390c0888a1bSJDBetteridge 
391c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
392c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
393c0888a1bSJDBetteridge 
394c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
395c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
396c0888a1bSJDBetteridge 
397c0888a1bSJDBetteridge   PetscCall(Sixth(&truth, 9));
398c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
399c0888a1bSJDBetteridge 
400c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
401c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
402c0888a1bSJDBetteridge   PetscCall(PetscFree(truth));
403c0888a1bSJDBetteridge 
404*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
405c0888a1bSJDBetteridge }
406c0888a1bSJDBetteridge 
407c0888a1bSJDBetteridge /* Tests functionality when setb is empty and setb is large */
408c0888a1bSJDBetteridge PetscErrorCode test_big_empty()
409c0888a1bSJDBetteridge {
410c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
411c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
412c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
413c0888a1bSJDBetteridge 
414c0888a1bSJDBetteridge   PetscFunctionBeginUser;
415c0888a1bSJDBetteridge 
416c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
417c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
418c0888a1bSJDBetteridge 
419c0888a1bSJDBetteridge   set_b[0] = 0;
420c0888a1bSJDBetteridge 
421c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
422c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
423c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
424c0888a1bSJDBetteridge 
425c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
426c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
427c0888a1bSJDBetteridge 
428*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429c0888a1bSJDBetteridge }
430c0888a1bSJDBetteridge 
431c0888a1bSJDBetteridge /* Tests functionality when setb is small and setb is large */
432c0888a1bSJDBetteridge PetscErrorCode test_big_small()
433c0888a1bSJDBetteridge {
434c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
435c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 8};
436c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
437c0888a1bSJDBetteridge 
438c0888a1bSJDBetteridge   PetscFunctionBeginUser;
439c0888a1bSJDBetteridge 
440c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
441c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_b));
442c0888a1bSJDBetteridge 
443c0888a1bSJDBetteridge   set_b[0] = 4;
444c0888a1bSJDBetteridge   set_b[1] = 1;
445c0888a1bSJDBetteridge   set_b[2] = 4;
446c0888a1bSJDBetteridge   set_b[3] = 8;
447c0888a1bSJDBetteridge   set_b[4] = 9;
448c0888a1bSJDBetteridge 
449c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
450c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
451c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
452c0888a1bSJDBetteridge 
453c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
454c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
455c0888a1bSJDBetteridge 
456*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
457c0888a1bSJDBetteridge }
458c0888a1bSJDBetteridge 
459c0888a1bSJDBetteridge /* Tests functionality when setb is medium sized and setb is large */
460c0888a1bSJDBetteridge PetscErrorCode test_big_moderate()
461c0888a1bSJDBetteridge {
462c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
463c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 8};
464c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
465c0888a1bSJDBetteridge 
466c0888a1bSJDBetteridge   PetscFunctionBeginUser;
467c0888a1bSJDBetteridge 
468c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
469c0888a1bSJDBetteridge   PetscCall(Fibonnaci(&set_b, 49));
470c0888a1bSJDBetteridge 
471c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
472c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
473c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
474c0888a1bSJDBetteridge 
475c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
476c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
477c0888a1bSJDBetteridge 
478*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
479c0888a1bSJDBetteridge }
480c0888a1bSJDBetteridge 
481c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large, in the opposite
482c0888a1bSJDBetteridge  order to test_big_big() */
483c0888a1bSJDBetteridge PetscErrorCode test_big_big_reversed()
484c0888a1bSJDBetteridge {
485c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
486c0888a1bSJDBetteridge   PetscInt64 *truth;
487c0888a1bSJDBetteridge   PetscMPIInt length = 1;
488c0888a1bSJDBetteridge 
489c0888a1bSJDBetteridge   PetscFunctionBeginUser;
490c0888a1bSJDBetteridge 
491c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
492c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
493c0888a1bSJDBetteridge 
494c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
495c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
496c0888a1bSJDBetteridge 
497c0888a1bSJDBetteridge   PetscCall(Sixth(&truth, 9));
498c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
499c0888a1bSJDBetteridge 
500c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
501c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
502c0888a1bSJDBetteridge   PetscCall(PetscFree(truth));
503c0888a1bSJDBetteridge 
504*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
505c0888a1bSJDBetteridge }
506c0888a1bSJDBetteridge 
507c0888a1bSJDBetteridge /* Main executes the individual tests in a predefined order */
508c0888a1bSJDBetteridge int main(int argc, char **argv)
509c0888a1bSJDBetteridge {
510c0888a1bSJDBetteridge   PetscFunctionBeginUser;
511c0888a1bSJDBetteridge   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
512c0888a1bSJDBetteridge 
513c0888a1bSJDBetteridge   /* Small tests */
514c0888a1bSJDBetteridge   /* Test different edge cases with small sets */
515c0888a1bSJDBetteridge   PetscCall(test_empty_empty());
516c0888a1bSJDBetteridge   PetscCall(test_a_empty());
517c0888a1bSJDBetteridge   PetscCall(test_b_empty());
518c0888a1bSJDBetteridge   PetscCall(test_identical());
519c0888a1bSJDBetteridge   PetscCall(test_disjoint());
520c0888a1bSJDBetteridge   PetscCall(test_single_common());
521c0888a1bSJDBetteridge   PetscCall(test_issue_1247());
522c0888a1bSJDBetteridge 
523c0888a1bSJDBetteridge   /* Big tests */
524c0888a1bSJDBetteridge   /* Test different edge cases with big sets */
525c0888a1bSJDBetteridge   PetscCall(test_empty_big());
526c0888a1bSJDBetteridge   PetscCall(test_small_big());
527c0888a1bSJDBetteridge   PetscCall(test_moderate_big());
528c0888a1bSJDBetteridge   PetscCall(test_big_big());
529c0888a1bSJDBetteridge   PetscCall(test_big_empty());
530c0888a1bSJDBetteridge   PetscCall(test_big_small());
531c0888a1bSJDBetteridge   PetscCall(test_big_moderate());
532c0888a1bSJDBetteridge   PetscCall(test_big_big_reversed());
533c0888a1bSJDBetteridge 
534c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ALL PASSED\n"));
535c0888a1bSJDBetteridge   PetscCall(PetscFinalize());
536c0888a1bSJDBetteridge   return 0;
537c0888a1bSJDBetteridge }
538c0888a1bSJDBetteridge 
539c0888a1bSJDBetteridge /*TEST
540c0888a1bSJDBetteridge 
541c0888a1bSJDBetteridge    test:
542758f5028SMatthew G. Knepley      suffix: 0
543c0888a1bSJDBetteridge 
544c0888a1bSJDBetteridge TEST*/
545