xref: /petsc/src/sys/tests/ex62.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
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 */
Fibonnaci(PetscInt64 ** set,PetscInt n)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;
3157508eceSPierre Jolivet   PetscAssert(n < 50, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "n must be less than 50");
32c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
33300f1712SStefano Zampini   (*set)[0] = n;
34*ac530a7eSPierre Jolivet   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 */
Square(PetscInt64 ** set,PetscInt n)39c0888a1bSJDBetteridge PetscErrorCode Square(PetscInt64 **set, PetscInt n)
40c0888a1bSJDBetteridge {
41c0888a1bSJDBetteridge   PetscInt64 ii;
42c0888a1bSJDBetteridge 
43c0888a1bSJDBetteridge   PetscFunctionBeginUser;
44c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
45300f1712SStefano Zampini   (*set)[0] = n;
46*ac530a7eSPierre Jolivet   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 */
Cube(PetscInt64 ** set,PetscInt n)51c0888a1bSJDBetteridge PetscErrorCode Cube(PetscInt64 **set, PetscInt n)
52c0888a1bSJDBetteridge {
53c0888a1bSJDBetteridge   PetscInt64 ii;
54c0888a1bSJDBetteridge 
55c0888a1bSJDBetteridge   PetscFunctionBeginUser;
56c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
57300f1712SStefano Zampini   (*set)[0] = n;
58*ac530a7eSPierre Jolivet   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 */
Sixth(PetscInt64 ** set,PetscInt n)63c0888a1bSJDBetteridge PetscErrorCode Sixth(PetscInt64 **set, PetscInt n)
64c0888a1bSJDBetteridge {
65c0888a1bSJDBetteridge   PetscInt64 ii;
66c0888a1bSJDBetteridge 
67c0888a1bSJDBetteridge   PetscFunctionBeginUser;
68c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(n + 1, set));
69300f1712SStefano Zampini   (*set)[0] = n;
70*ac530a7eSPierre Jolivet   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 */
PrintSet(PetscInt64 * set)75c0888a1bSJDBetteridge PetscErrorCode PrintSet(PetscInt64 *set)
76c0888a1bSJDBetteridge {
77c0888a1bSJDBetteridge   char text[64];
78c0888a1bSJDBetteridge 
79c0888a1bSJDBetteridge   PetscFunctionBeginUser;
80c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "["));
81300f1712SStefano Zampini   for (PetscInt64 ii = 1; ii <= set[0]; ii++) {
82c0888a1bSJDBetteridge     PetscCall(PetscFormatConvert(" %" PetscInt64_FMT ",", text));
83c0888a1bSJDBetteridge     PetscCall(PetscPrintf(PETSC_COMM_WORLD, text, set[ii]));
84c0888a1bSJDBetteridge   }
85c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n"));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c0888a1bSJDBetteridge }
88c0888a1bSJDBetteridge 
89c0888a1bSJDBetteridge /* Check set equality */
AssertSetsEqual(PetscInt64 * set,PetscInt64 * true_set)90c0888a1bSJDBetteridge PetscErrorCode AssertSetsEqual(PetscInt64 *set, PetscInt64 *true_set)
91c0888a1bSJDBetteridge {
92c0888a1bSJDBetteridge   PetscInt ii;
93c0888a1bSJDBetteridge 
94c0888a1bSJDBetteridge   PetscFunctionBeginUser;
9557508eceSPierre Jolivet   PetscAssert(set[0] == true_set[0], PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets of different sizes");
9657508eceSPierre 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");
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98c0888a1bSJDBetteridge }
99c0888a1bSJDBetteridge 
100c0888a1bSJDBetteridge /* Tests functionality when two enpty sets are passed */
test_empty_empty()101c0888a1bSJDBetteridge PetscErrorCode test_empty_empty()
102c0888a1bSJDBetteridge {
103c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
104c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
105c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
106c0888a1bSJDBetteridge 
107c0888a1bSJDBetteridge   PetscFunctionBeginUser;
108c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
109c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
110c0888a1bSJDBetteridge 
111c0888a1bSJDBetteridge   set_a[0] = 0;
112c0888a1bSJDBetteridge 
113c0888a1bSJDBetteridge   set_b[0] = 0;
114c0888a1bSJDBetteridge 
115c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
116c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
117c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
118c0888a1bSJDBetteridge 
119c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
120c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122c0888a1bSJDBetteridge }
123c0888a1bSJDBetteridge 
124c0888a1bSJDBetteridge /* Tests functionality when seta is empty */
test_a_empty()125c0888a1bSJDBetteridge PetscErrorCode test_a_empty()
126c0888a1bSJDBetteridge {
127c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
128c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
129c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
130c0888a1bSJDBetteridge 
131c0888a1bSJDBetteridge   PetscFunctionBeginUser;
132c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
133c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_b));
134c0888a1bSJDBetteridge 
135c0888a1bSJDBetteridge   set_a[0] = 0;
136c0888a1bSJDBetteridge 
137c0888a1bSJDBetteridge   set_b[0] = 1;
138c0888a1bSJDBetteridge   set_b[1] = 1;
139c0888a1bSJDBetteridge 
140c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
141c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
142c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
143c0888a1bSJDBetteridge 
144c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
145c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147c0888a1bSJDBetteridge }
148c0888a1bSJDBetteridge 
149c0888a1bSJDBetteridge /* Tests functionality when setb is empty */
test_b_empty()150c0888a1bSJDBetteridge PetscErrorCode test_b_empty()
151c0888a1bSJDBetteridge {
152c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
153c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
154c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
155c0888a1bSJDBetteridge 
156c0888a1bSJDBetteridge   PetscFunctionBeginUser;
157c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_a));
158c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
159c0888a1bSJDBetteridge 
160c0888a1bSJDBetteridge   set_a[0] = 1;
161c0888a1bSJDBetteridge   set_a[1] = 1;
162c0888a1bSJDBetteridge 
163c0888a1bSJDBetteridge   set_b[0] = 0;
164c0888a1bSJDBetteridge 
165c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
166c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
167c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
168c0888a1bSJDBetteridge 
169c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
170c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c0888a1bSJDBetteridge }
173c0888a1bSJDBetteridge 
174c0888a1bSJDBetteridge /* Tests functionality when both sets are identical */
test_identical()175c0888a1bSJDBetteridge PetscErrorCode test_identical()
176c0888a1bSJDBetteridge {
177c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
178c0888a1bSJDBetteridge   PetscInt64  truth[] = {3, 1, 4, 9};
179c0888a1bSJDBetteridge   PetscMPIInt length  = 4;
180c0888a1bSJDBetteridge 
181c0888a1bSJDBetteridge   PetscFunctionBeginUser;
182c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
183c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_b));
184c0888a1bSJDBetteridge 
185c0888a1bSJDBetteridge   set_a[0] = 3;
186c0888a1bSJDBetteridge   set_a[1] = 1;
187c0888a1bSJDBetteridge   set_a[2] = 4;
188c0888a1bSJDBetteridge   set_a[3] = 9;
189c0888a1bSJDBetteridge 
190c0888a1bSJDBetteridge   set_b[0] = 3;
191c0888a1bSJDBetteridge   set_b[1] = 1;
192c0888a1bSJDBetteridge   set_b[2] = 4;
193c0888a1bSJDBetteridge   set_b[3] = 9;
194c0888a1bSJDBetteridge 
195c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
196c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
197c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
198c0888a1bSJDBetteridge 
199c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
200c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202c0888a1bSJDBetteridge }
203c0888a1bSJDBetteridge 
204c0888a1bSJDBetteridge /* Tests functionality when sets have no elements in common */
test_disjoint()205c0888a1bSJDBetteridge PetscErrorCode test_disjoint()
206c0888a1bSJDBetteridge {
207c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
208c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
209c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
210c0888a1bSJDBetteridge 
211c0888a1bSJDBetteridge   PetscFunctionBeginUser;
212c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
213c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_b));
214c0888a1bSJDBetteridge 
215c0888a1bSJDBetteridge   set_a[0] = 3;
216c0888a1bSJDBetteridge   set_a[1] = 1;
217c0888a1bSJDBetteridge   set_a[2] = 4;
218c0888a1bSJDBetteridge   set_a[3] = 9;
219c0888a1bSJDBetteridge 
220c0888a1bSJDBetteridge   set_b[0] = 3;
221c0888a1bSJDBetteridge   set_b[1] = 2;
222c0888a1bSJDBetteridge   set_b[2] = 6;
223c0888a1bSJDBetteridge   set_b[3] = 8;
224c0888a1bSJDBetteridge 
225c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
226c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
227c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
228c0888a1bSJDBetteridge 
229c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
230c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232c0888a1bSJDBetteridge }
233c0888a1bSJDBetteridge 
234c0888a1bSJDBetteridge /* Tests functionality when sets only have one element in common */
test_single_common()235c0888a1bSJDBetteridge PetscErrorCode test_single_common()
236c0888a1bSJDBetteridge {
237c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
238c0888a1bSJDBetteridge   PetscInt64  truth[] = {1, 4};
239c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
240c0888a1bSJDBetteridge 
241c0888a1bSJDBetteridge   PetscFunctionBeginUser;
242c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(4, &set_a));
243c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_b));
244c0888a1bSJDBetteridge 
245c0888a1bSJDBetteridge   set_a[0] = 3;
246c0888a1bSJDBetteridge   set_a[1] = 1;
247c0888a1bSJDBetteridge   set_a[2] = 4;
248c0888a1bSJDBetteridge   set_a[3] = 9;
249c0888a1bSJDBetteridge 
250c0888a1bSJDBetteridge   set_b[0] = 3;
251c0888a1bSJDBetteridge   set_b[1] = 2;
252c0888a1bSJDBetteridge   set_b[2] = 4;
253c0888a1bSJDBetteridge   set_b[3] = 6;
254c0888a1bSJDBetteridge   set_b[4] = 8;
255c0888a1bSJDBetteridge 
256c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
257c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
258c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
259c0888a1bSJDBetteridge 
260c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
261c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263c0888a1bSJDBetteridge }
264c0888a1bSJDBetteridge 
265c0888a1bSJDBetteridge /* Specific test case flagged by PETSc issue #1247 */
test_issue_1247()266c0888a1bSJDBetteridge PetscErrorCode test_issue_1247()
267c0888a1bSJDBetteridge {
268c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
269c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
270c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
271c0888a1bSJDBetteridge 
272c0888a1bSJDBetteridge   PetscFunctionBeginUser;
273c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(3, &set_a));
274c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(2, &set_b));
275c0888a1bSJDBetteridge 
276c0888a1bSJDBetteridge   set_a[0] = 2;
277c0888a1bSJDBetteridge   set_a[1] = 2;
278c0888a1bSJDBetteridge   set_a[2] = 3;
279c0888a1bSJDBetteridge 
280c0888a1bSJDBetteridge   set_b[0] = 1;
281c0888a1bSJDBetteridge   set_b[1] = 1;
282c0888a1bSJDBetteridge 
283c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
284c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
285c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
286c0888a1bSJDBetteridge 
287c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
288c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
290c0888a1bSJDBetteridge }
291c0888a1bSJDBetteridge 
292c0888a1bSJDBetteridge /* Tests functionality when seta is empty and setb is large */
test_empty_big()293c0888a1bSJDBetteridge PetscErrorCode test_empty_big()
294c0888a1bSJDBetteridge {
295c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
296c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
297c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
298c0888a1bSJDBetteridge 
299c0888a1bSJDBetteridge   PetscFunctionBeginUser;
300c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_a));
301c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
302c0888a1bSJDBetteridge 
303c0888a1bSJDBetteridge   set_a[0] = 0;
304c0888a1bSJDBetteridge 
305c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
306c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
307c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
308c0888a1bSJDBetteridge 
309c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
310c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312c0888a1bSJDBetteridge }
313c0888a1bSJDBetteridge 
314c0888a1bSJDBetteridge /* Tests functionality when seta is small and setb is large */
test_small_big()315c0888a1bSJDBetteridge PetscErrorCode test_small_big()
316c0888a1bSJDBetteridge {
317c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
318c0888a1bSJDBetteridge   PetscInt64  truth[] = {3, 1, 4, 9};
319c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
320c0888a1bSJDBetteridge 
321c0888a1bSJDBetteridge   PetscFunctionBeginUser;
322c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_a));
323c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
324c0888a1bSJDBetteridge 
325c0888a1bSJDBetteridge   set_a[0] = 4;
326c0888a1bSJDBetteridge   set_a[1] = 1;
327c0888a1bSJDBetteridge   set_a[2] = 4;
328c0888a1bSJDBetteridge   set_a[3] = 8;
329c0888a1bSJDBetteridge   set_a[4] = 9;
330c0888a1bSJDBetteridge 
331c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
332c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
333c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
334c0888a1bSJDBetteridge 
335c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
336c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338c0888a1bSJDBetteridge }
339c0888a1bSJDBetteridge 
340c0888a1bSJDBetteridge /* Tests functionality when seta is medium sized and setb is large */
test_moderate_big()341c0888a1bSJDBetteridge PetscErrorCode test_moderate_big()
342c0888a1bSJDBetteridge {
343c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
344c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 144};
345c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
346c0888a1bSJDBetteridge 
347c0888a1bSJDBetteridge   PetscFunctionBeginUser;
348c0888a1bSJDBetteridge   PetscCall(Fibonnaci(&set_a, 49));
349c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
350c0888a1bSJDBetteridge 
351c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
352c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
353c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
354c0888a1bSJDBetteridge 
355c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
356c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358c0888a1bSJDBetteridge }
359c0888a1bSJDBetteridge 
360c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large */
test_big_big()361c0888a1bSJDBetteridge PetscErrorCode test_big_big()
362c0888a1bSJDBetteridge {
363c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
364c0888a1bSJDBetteridge   PetscInt64 *truth;
365c0888a1bSJDBetteridge   PetscMPIInt length = 1;
366c0888a1bSJDBetteridge 
367c0888a1bSJDBetteridge   PetscFunctionBeginUser;
368c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
369c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
370c0888a1bSJDBetteridge 
371c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
372c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
373c0888a1bSJDBetteridge 
374c0888a1bSJDBetteridge   PetscCall(Sixth(&truth, 9));
375c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
376c0888a1bSJDBetteridge 
377c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
378c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
379c0888a1bSJDBetteridge   PetscCall(PetscFree(truth));
3803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
381c0888a1bSJDBetteridge }
382c0888a1bSJDBetteridge 
383c0888a1bSJDBetteridge /* Tests functionality when setb is empty and setb is large */
test_big_empty()384c0888a1bSJDBetteridge PetscErrorCode test_big_empty()
385c0888a1bSJDBetteridge {
386c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
387c0888a1bSJDBetteridge   PetscInt64  truth[] = {0};
388c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
389c0888a1bSJDBetteridge 
390c0888a1bSJDBetteridge   PetscFunctionBeginUser;
391c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
392c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(1, &set_b));
393c0888a1bSJDBetteridge 
394c0888a1bSJDBetteridge   set_b[0] = 0;
395c0888a1bSJDBetteridge 
396c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
397c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
398c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
399c0888a1bSJDBetteridge 
400c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
401c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
403c0888a1bSJDBetteridge }
404c0888a1bSJDBetteridge 
405c0888a1bSJDBetteridge /* Tests functionality when setb is small and setb is large */
test_big_small()406c0888a1bSJDBetteridge PetscErrorCode test_big_small()
407c0888a1bSJDBetteridge {
408c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
409c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 8};
410c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
411c0888a1bSJDBetteridge 
412c0888a1bSJDBetteridge   PetscFunctionBeginUser;
413c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
414c0888a1bSJDBetteridge   PetscCall(PetscMalloc1(5, &set_b));
415c0888a1bSJDBetteridge 
416c0888a1bSJDBetteridge   set_b[0] = 4;
417c0888a1bSJDBetteridge   set_b[1] = 1;
418c0888a1bSJDBetteridge   set_b[2] = 4;
419c0888a1bSJDBetteridge   set_b[3] = 8;
420c0888a1bSJDBetteridge   set_b[4] = 9;
421c0888a1bSJDBetteridge 
422c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
423c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
424c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
425c0888a1bSJDBetteridge 
426c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
427c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429c0888a1bSJDBetteridge }
430c0888a1bSJDBetteridge 
431c0888a1bSJDBetteridge /* Tests functionality when setb is medium sized and setb is large */
test_big_moderate()432c0888a1bSJDBetteridge PetscErrorCode test_big_moderate()
433c0888a1bSJDBetteridge {
434c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
435c0888a1bSJDBetteridge   PetscInt64  truth[] = {2, 1, 8};
436c0888a1bSJDBetteridge   PetscMPIInt length  = 1;
437c0888a1bSJDBetteridge 
438c0888a1bSJDBetteridge   PetscFunctionBeginUser;
439c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
440c0888a1bSJDBetteridge   PetscCall(Fibonnaci(&set_b, 49));
441c0888a1bSJDBetteridge 
442c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
443c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
444c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
445c0888a1bSJDBetteridge 
446c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
447c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
449c0888a1bSJDBetteridge }
450c0888a1bSJDBetteridge 
451c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large, in the opposite
452c0888a1bSJDBetteridge  order to test_big_big() */
test_big_big_reversed()453c0888a1bSJDBetteridge PetscErrorCode test_big_big_reversed()
454c0888a1bSJDBetteridge {
455c0888a1bSJDBetteridge   PetscInt64 *set_a, *set_b;
456c0888a1bSJDBetteridge   PetscInt64 *truth;
457c0888a1bSJDBetteridge   PetscMPIInt length = 1;
458c0888a1bSJDBetteridge 
459c0888a1bSJDBetteridge   PetscFunctionBeginUser;
460c0888a1bSJDBetteridge   PetscCall(Cube(&set_a, 999));
461c0888a1bSJDBetteridge   PetscCall(Square(&set_b, 999));
462c0888a1bSJDBetteridge 
463c0888a1bSJDBetteridge   PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL);
464c0888a1bSJDBetteridge   PetscCall(PrintSet(set_a));
465c0888a1bSJDBetteridge 
466c0888a1bSJDBetteridge   PetscCall(Sixth(&truth, 9));
467c0888a1bSJDBetteridge   PetscCall(AssertSetsEqual(set_a, truth));
468c0888a1bSJDBetteridge 
469c0888a1bSJDBetteridge   PetscCall(PetscFree(set_a));
470c0888a1bSJDBetteridge   PetscCall(PetscFree(set_b));
471c0888a1bSJDBetteridge   PetscCall(PetscFree(truth));
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473c0888a1bSJDBetteridge }
474c0888a1bSJDBetteridge 
475c0888a1bSJDBetteridge /* Main executes the individual tests in a predefined order */
main(int argc,char ** argv)476c0888a1bSJDBetteridge int main(int argc, char **argv)
477c0888a1bSJDBetteridge {
478c0888a1bSJDBetteridge   PetscFunctionBeginUser;
479c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
480c0888a1bSJDBetteridge 
481c0888a1bSJDBetteridge   /* Small tests */
482c0888a1bSJDBetteridge   /* Test different edge cases with small sets */
483c0888a1bSJDBetteridge   PetscCall(test_empty_empty());
484c0888a1bSJDBetteridge   PetscCall(test_a_empty());
485c0888a1bSJDBetteridge   PetscCall(test_b_empty());
486c0888a1bSJDBetteridge   PetscCall(test_identical());
487c0888a1bSJDBetteridge   PetscCall(test_disjoint());
488c0888a1bSJDBetteridge   PetscCall(test_single_common());
489c0888a1bSJDBetteridge   PetscCall(test_issue_1247());
490c0888a1bSJDBetteridge 
491c0888a1bSJDBetteridge   /* Big tests */
492c0888a1bSJDBetteridge   /* Test different edge cases with big sets */
493c0888a1bSJDBetteridge   PetscCall(test_empty_big());
494c0888a1bSJDBetteridge   PetscCall(test_small_big());
495c0888a1bSJDBetteridge   PetscCall(test_moderate_big());
496c0888a1bSJDBetteridge   PetscCall(test_big_big());
497c0888a1bSJDBetteridge   PetscCall(test_big_empty());
498c0888a1bSJDBetteridge   PetscCall(test_big_small());
499c0888a1bSJDBetteridge   PetscCall(test_big_moderate());
500c0888a1bSJDBetteridge   PetscCall(test_big_big_reversed());
501c0888a1bSJDBetteridge 
502c0888a1bSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ALL PASSED\n"));
503c0888a1bSJDBetteridge   PetscCall(PetscFinalize());
504c0888a1bSJDBetteridge   return 0;
505c0888a1bSJDBetteridge }
506c0888a1bSJDBetteridge 
507c0888a1bSJDBetteridge /*TEST
508c0888a1bSJDBetteridge 
509c0888a1bSJDBetteridge    test:
510758f5028SMatthew G. Knepley      suffix: 0
511c0888a1bSJDBetteridge 
512c0888a1bSJDBetteridge TEST*/
513