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