xref: /petsc/src/sys/utils/sorti.c (revision 0619917b5a674bb687c64e7daba2ab22be99af31)
1 
2 /*
3    This file contains routines for sorting integers. Values are sorted in place.
4    One can use src/sys/tests/ex52.c for benchmarking.
5  */
6 #include <petsc/private/petscimpl.h> /*I  "petscsys.h"  I*/
7 #include <petsc/private/hashseti.h>
8 
9 #define MEDIAN3(v, a, b, c) (v[a] < v[b] ? (v[b] < v[c] ? (b) : (v[a] < v[c] ? (c) : (a))) : (v[c] < v[b] ? (b) : (v[a] < v[c] ? (a) : (c))))
10 
11 #define MEDIAN(v, right) MEDIAN3(v, right / 4, right / 2, right / 4 * 3)
12 
13 /* Swap one, two or three pairs. Each pair can have its own type */
14 #define SWAP1(a, b, t1) \
15   do { \
16     t1 = a; \
17     a  = b; \
18     b  = t1; \
19   } while (0)
20 #define SWAP2(a, b, c, d, t1, t2) \
21   do { \
22     t1 = a; \
23     a  = b; \
24     b  = t1; \
25     t2 = c; \
26     c  = d; \
27     d  = t2; \
28   } while (0)
29 #define SWAP3(a, b, c, d, e, f, t1, t2, t3) \
30   do { \
31     t1 = a; \
32     a  = b; \
33     b  = t1; \
34     t2 = c; \
35     c  = d; \
36     d  = t2; \
37     t3 = e; \
38     e  = f; \
39     f  = t3; \
40   } while (0)
41 
42 /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
43 #define SWAP2Data(a, b, c, d, t1, t2, siz) \
44   do { \
45     t1 = a; \
46     a  = b; \
47     b  = t1; \
48     PetscCall(PetscMemcpy(t2, c, siz)); \
49     PetscCall(PetscMemcpy(c, d, siz)); \
50     PetscCall(PetscMemcpy(d, t2, siz)); \
51   } while (0)
52 
53 /*
54    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
55 
56    Input Parameters:
57     + X         - array to partition
58     . pivot     - a pivot of X[]
59     . t1        - temp variable for X
60     - lo,hi     - lower and upper bound of the array
61 
62    Output Parameters:
63     + l,r       - of type PetscInt
64 
65    Note:
66     The TwoWayPartition2/3 variants also partition other arrays along with X.
67     These arrays can have different types, so they provide their own temp t2,t3
68  */
69 #define TwoWayPartition1(X, pivot, t1, lo, hi, l, r) \
70   do { \
71     l = lo; \
72     r = hi; \
73     while (1) { \
74       while (X[l] < pivot) l++; \
75       while (X[r] > pivot) r--; \
76       if (l >= r) { \
77         r++; \
78         break; \
79       } \
80       SWAP1(X[l], X[r], t1); \
81       l++; \
82       r--; \
83     } \
84   } while (0)
85 
86 /*
87    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
88 
89    Input Parameters:
90     + X         - array to partition
91     . pivot     - a pivot of X[]
92     . t1        - temp variable for X
93     - lo,hi     - lower and upper bound of the array
94 
95    Output Parameters:
96     + l,r       - of type PetscInt
97 
98    Note:
99     The TwoWayPartition2/3 variants also partition other arrays along with X.
100     These arrays can have different types, so they provide their own temp t2,t3
101  */
102 #define TwoWayPartitionReverse1(X, pivot, t1, lo, hi, l, r) \
103   do { \
104     l = lo; \
105     r = hi; \
106     while (1) { \
107       while (X[l] > pivot) l++; \
108       while (X[r] < pivot) r--; \
109       if (l >= r) { \
110         r++; \
111         break; \
112       } \
113       SWAP1(X[l], X[r], t1); \
114       l++; \
115       r--; \
116     } \
117   } while (0)
118 
119 #define TwoWayPartition2(X, Y, pivot, t1, t2, lo, hi, l, r) \
120   do { \
121     l = lo; \
122     r = hi; \
123     while (1) { \
124       while (X[l] < pivot) l++; \
125       while (X[r] > pivot) r--; \
126       if (l >= r) { \
127         r++; \
128         break; \
129       } \
130       SWAP2(X[l], X[r], Y[l], Y[r], t1, t2); \
131       l++; \
132       r--; \
133     } \
134   } while (0)
135 
136 #define TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, lo, hi, l, r) \
137   do { \
138     l = lo; \
139     r = hi; \
140     while (1) { \
141       while (X[l] < pivot) l++; \
142       while (X[r] > pivot) r--; \
143       if (l >= r) { \
144         r++; \
145         break; \
146       } \
147       SWAP3(X[l], X[r], Y[l], Y[r], Z[l], Z[r], t1, t2, t3); \
148       l++; \
149       r--; \
150     } \
151   } while (0)
152 
153 /* Templates for similar functions used below */
154 #define QuickSort1(FuncName, X, n, pivot, t1) \
155   do { \
156     PetscCount i, j, p, l, r, hi = n - 1; \
157     if (n < 8) { \
158       for (i = 0; i < n; i++) { \
159         pivot = X[i]; \
160         for (j = i + 1; j < n; j++) { \
161           if (pivot > X[j]) { \
162             SWAP1(X[i], X[j], t1); \
163             pivot = X[i]; \
164           } \
165         } \
166       } \
167     } else { \
168       p     = MEDIAN(X, hi); \
169       pivot = X[p]; \
170       TwoWayPartition1(X, pivot, t1, 0, hi, l, r); \
171       PetscCall(FuncName(l, X)); \
172       PetscCall(FuncName(hi - r + 1, X + r)); \
173     } \
174   } while (0)
175 
176 /* Templates for similar functions used below */
177 #define QuickSortReverse1(FuncName, X, n, pivot, t1) \
178   do { \
179     PetscCount i, j, p, l, r, hi = n - 1; \
180     if (n < 8) { \
181       for (i = 0; i < n; i++) { \
182         pivot = X[i]; \
183         for (j = i + 1; j < n; j++) { \
184           if (pivot < X[j]) { \
185             SWAP1(X[i], X[j], t1); \
186             pivot = X[i]; \
187           } \
188         } \
189       } \
190     } else { \
191       p     = MEDIAN(X, hi); \
192       pivot = X[p]; \
193       TwoWayPartitionReverse1(X, pivot, t1, 0, hi, l, r); \
194       PetscCall(FuncName(l, X)); \
195       PetscCall(FuncName(hi - r + 1, X + r)); \
196     } \
197   } while (0)
198 
199 #define QuickSort2(FuncName, X, Y, n, pivot, t1, t2) \
200   do { \
201     PetscCount i, j, p, l, r, hi = n - 1; \
202     if (n < 8) { \
203       for (i = 0; i < n; i++) { \
204         pivot = X[i]; \
205         for (j = i + 1; j < n; j++) { \
206           if (pivot > X[j]) { \
207             SWAP2(X[i], X[j], Y[i], Y[j], t1, t2); \
208             pivot = X[i]; \
209           } \
210         } \
211       } \
212     } else { \
213       p     = MEDIAN(X, hi); \
214       pivot = X[p]; \
215       TwoWayPartition2(X, Y, pivot, t1, t2, 0, hi, l, r); \
216       PetscCall(FuncName(l, X, Y)); \
217       PetscCall(FuncName(hi - r + 1, X + r, Y + r)); \
218     } \
219   } while (0)
220 
221 #define QuickSort3(FuncName, X, Y, Z, n, pivot, t1, t2, t3) \
222   do { \
223     PetscCount i, j, p, l, r, hi = n - 1; \
224     if (n < 8) { \
225       for (i = 0; i < n; i++) { \
226         pivot = X[i]; \
227         for (j = i + 1; j < n; j++) { \
228           if (pivot > X[j]) { \
229             SWAP3(X[i], X[j], Y[i], Y[j], Z[i], Z[j], t1, t2, t3); \
230             pivot = X[i]; \
231           } \
232         } \
233       } \
234     } else { \
235       p     = MEDIAN(X, hi); \
236       pivot = X[p]; \
237       TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, 0, hi, l, r); \
238       PetscCall(FuncName(l, X, Y, Z)); \
239       PetscCall(FuncName(hi - r + 1, X + r, Y + r, Z + r)); \
240     } \
241   } while (0)
242 
243 /*@
244   PetscSortedInt - Determines whether the `PetscInt` array is sorted.
245 
246   Not Collective
247 
248   Input Parameters:
249 + n - number of values
250 - X - array of integers
251 
252   Output Parameter:
253 . sorted - flag whether the array is sorted
254 
255   Level: intermediate
256 
257 .seealso: `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
258 @*/
259 PetscErrorCode PetscSortedInt(PetscInt n, const PetscInt X[], PetscBool *sorted)
260 {
261   PetscFunctionBegin;
262   if (n) PetscAssertPointer(X, 2);
263   PetscAssertPointer(sorted, 3);
264   PetscSorted(n, X, *sorted);
265   PetscFunctionReturn(PETSC_SUCCESS);
266 }
267 
268 /*@
269   PetscSortedInt64 - Determines whether the `PetscInt64` array is sorted.
270 
271   Not Collective
272 
273   Input Parameters:
274 + n - number of values
275 - X - array of integers
276 
277   Output Parameter:
278 . sorted - flag whether the array is sorted
279 
280   Level: intermediate
281 
282 .seealso: `PetscSortInt64()`, `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
283 @*/
284 PetscErrorCode PetscSortedInt64(PetscInt n, const PetscInt64 X[], PetscBool *sorted)
285 {
286   PetscFunctionBegin;
287   if (n) PetscAssertPointer(X, 2);
288   PetscAssertPointer(sorted, 3);
289   PetscSorted(n, X, *sorted);
290   PetscFunctionReturn(PETSC_SUCCESS);
291 }
292 
293 /*@
294   PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.
295 
296   Not Collective
297 
298   Input Parameters:
299 + n - number of values
300 - X - array of integers
301 
302   Note:
303   This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
304   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
305   code to see which routine is fastest.
306 
307   Level: intermediate
308 
309 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
310 @*/
311 PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[])
312 {
313   PetscInt pivot, t1;
314 
315   PetscFunctionBegin;
316   if (n) PetscAssertPointer(X, 2);
317   QuickSort1(PetscSortInt, X, n, pivot, t1);
318   PetscFunctionReturn(PETSC_SUCCESS);
319 }
320 
321 /*@
322   PetscSortInt64 - Sorts an array of `PetscInt64` in place in increasing order.
323 
324   Not Collective
325 
326   Input Parameters:
327 + n - number of values
328 - X - array of integers
329 
330   Notes:
331   This function sorts `PetscCount`s assumed to be in completely random order
332 
333   Level: intermediate
334 
335 .seealso: `PetscSortInt()`
336 @*/
337 PetscErrorCode PetscSortInt64(PetscInt n, PetscInt64 X[])
338 {
339   PetscCount pivot, t1;
340 
341   PetscFunctionBegin;
342   if (n) PetscAssertPointer(X, 2);
343   QuickSort1(PetscSortInt64, X, n, pivot, t1);
344   PetscFunctionReturn(PETSC_SUCCESS);
345 }
346 
347 /*@
348   PetscSortCount - Sorts an array of integers in place in increasing order.
349 
350   Not Collective
351 
352   Input Parameters:
353 + n - number of values
354 - X - array of integers
355 
356   Notes:
357   This function sorts `PetscCount`s assumed to be in completely random order
358 
359   Level: intermediate
360 
361 .seealso: `PetscSortInt()`
362 @*/
363 PetscErrorCode PetscSortCount(PetscInt n, PetscCount X[])
364 {
365   PetscCount pivot, t1;
366 
367   PetscFunctionBegin;
368   if (n) PetscAssertPointer(X, 2);
369   QuickSort1(PetscSortCount, X, n, pivot, t1);
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372 
373 /*@
374   PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
375 
376   Not Collective
377 
378   Input Parameters:
379 + n - number of values
380 - X - array of integers
381 
382   Level: intermediate
383 
384 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
385 @*/
386 PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[])
387 {
388   PetscInt pivot, t1;
389 
390   PetscFunctionBegin;
391   if (n) PetscAssertPointer(X, 2);
392   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
393   PetscFunctionReturn(PETSC_SUCCESS);
394 }
395 
396 /*@
397   PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array
398 
399   Not Collective
400 
401   Input Parameters:
402 + n - number of values
403 - X - sorted array of integers
404 
405   Output Parameter:
406 . n - number of non-redundant values
407 
408   Level: intermediate
409 
410 .seealso: `PetscSortInt()`
411 @*/
412 PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[])
413 {
414   PetscInt i, s = 0, N = *n, b = 0;
415 
416   PetscFunctionBegin;
417   PetscAssertPointer(n, 1);
418   PetscCheckSorted(*n, X);
419   for (i = 0; i < N - 1; i++) {
420     if (X[b + s + 1] != X[b]) {
421       X[b + 1] = X[b + s + 1];
422       b++;
423     } else s++;
424   }
425   *n = N - s;
426   PetscFunctionReturn(PETSC_SUCCESS);
427 }
428 
429 /*@
430   PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates
431 
432   Not Collective
433 
434   Input Parameters:
435 + n - number of values
436 - X - sorted array of integers
437 
438   Output Parameter:
439 . flg - True if the array has dups, otherwise false
440 
441   Level: intermediate
442 
443 .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
444 @*/
445 PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg)
446 {
447   PetscInt i;
448 
449   PetscFunctionBegin;
450   PetscCheckSorted(n, X);
451   *flg = PETSC_FALSE;
452   for (i = 0; i < n - 1; i++) {
453     if (X[i + 1] == X[i]) {
454       *flg = PETSC_TRUE;
455       break;
456     }
457   }
458   PetscFunctionReturn(PETSC_SUCCESS);
459 }
460 
461 /*@
462   PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
463 
464   Not Collective
465 
466   Input Parameters:
467 + n - number of values
468 - X - array of integers
469 
470   Output Parameter:
471 . n - number of non-redundant values
472 
473   Level: intermediate
474 
475 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
476 @*/
477 PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
478 {
479   PetscFunctionBegin;
480   PetscAssertPointer(n, 1);
481   PetscCall(PetscSortInt(*n, X));
482   PetscCall(PetscSortedRemoveDupsInt(n, X));
483   PetscFunctionReturn(PETSC_SUCCESS);
484 }
485 
486 /*@
487   PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt`
488 
489   Not Collective
490 
491   Input Parameters:
492 + key - the integer to locate
493 . n   - number of values in the array
494 - X   - array of integers
495 
496   Output Parameter:
497 . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
498 
499   Level: intermediate
500 
501 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
502 @*/
503 PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
504 {
505   PetscInt lo = 0, hi = n;
506 
507   PetscFunctionBegin;
508   PetscAssertPointer(loc, 4);
509   if (!n) {
510     *loc = -1;
511     PetscFunctionReturn(PETSC_SUCCESS);
512   }
513   PetscAssertPointer(X, 3);
514   PetscCheckSorted(n, X);
515   while (hi - lo > 1) {
516     PetscInt mid = lo + (hi - lo) / 2;
517     if (key < X[mid]) hi = mid;
518     else lo = mid;
519   }
520   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
521   PetscFunctionReturn(PETSC_SUCCESS);
522 }
523 
524 /*@
525   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
526 
527   Not Collective
528 
529   Input Parameters:
530 + n - number of values in the array
531 - X - array of integers
532 
533   Output Parameter:
534 . dups - True if the array has dups, otherwise false
535 
536   Level: intermediate
537 
538 .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
539 @*/
540 PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
541 {
542   PetscInt   i;
543   PetscHSetI ht;
544   PetscBool  missing;
545 
546   PetscFunctionBegin;
547   if (n) PetscAssertPointer(X, 2);
548   PetscAssertPointer(dups, 3);
549   *dups = PETSC_FALSE;
550   if (n > 1) {
551     PetscCall(PetscHSetICreate(&ht));
552     PetscCall(PetscHSetIResize(ht, n));
553     for (i = 0; i < n; i++) {
554       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
555       if (!missing) {
556         *dups = PETSC_TRUE;
557         break;
558       }
559     }
560     PetscCall(PetscHSetIDestroy(&ht));
561   }
562   PetscFunctionReturn(PETSC_SUCCESS);
563 }
564 
565 /*@
566   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
567 
568   Not Collective
569 
570   Input Parameters:
571 + key - the integer to locate
572 . n   - number of values in the array
573 - X   - array of integers
574 
575   Output Parameter:
576 . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
577 
578   Level: intermediate
579 
580 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
581 @*/
582 PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
583 {
584   PetscInt lo = 0, hi = n;
585 
586   PetscFunctionBegin;
587   PetscAssertPointer(loc, 4);
588   if (!n) {
589     *loc = -1;
590     PetscFunctionReturn(PETSC_SUCCESS);
591   }
592   PetscAssertPointer(X, 3);
593   PetscCheckSorted(n, X);
594   while (hi - lo > 1) {
595     PetscInt mid = lo + (hi - lo) / 2;
596     if (key < X[mid]) hi = mid;
597     else lo = mid;
598   }
599   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
600   PetscFunctionReturn(PETSC_SUCCESS);
601 }
602 
603 /*@
604   PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
605   changes a second array of `PetscInt` to match the sorted first array.
606 
607   Not Collective
608 
609   Input Parameters:
610 + n - number of values
611 . X - array of integers
612 - Y - second array of integers
613 
614   Level: intermediate
615 
616 .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
617 @*/
618 PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[])
619 {
620   PetscInt pivot, t1, t2;
621 
622   PetscFunctionBegin;
623   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
624   PetscFunctionReturn(PETSC_SUCCESS);
625 }
626 
627 /*@
628   PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
629   changes a pair of `PetscInt` arrays to match the sorted first array.
630 
631   Not Collective
632 
633   Input Parameters:
634 + n - number of values
635 . X - array of integers
636 . Y - second array of integers (first array of the pair)
637 - Z - third array of integers  (second array of the pair)
638 
639   Level: intermediate
640 
641 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
642 @*/
643 PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[])
644 {
645   PetscInt pivot, t1, t2, t3;
646 
647   PetscFunctionBegin;
648   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
649   PetscFunctionReturn(PETSC_SUCCESS);
650 }
651 
652 /*@
653   PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
654   changes a second array of `PetscCount` to match the sorted first array.
655 
656   Not Collective
657 
658   Input Parameters:
659 + n - number of values
660 . X - array of integers
661 - Y - second array of PetscCounts (signed integers)
662 
663   Level: intermediate
664 
665 .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
666 @*/
667 PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
668 {
669   PetscInt   pivot, t1;
670   PetscCount t2;
671 
672   PetscFunctionBegin;
673   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
674   PetscFunctionReturn(PETSC_SUCCESS);
675 }
676 
677 /*@
678   PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
679   changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
680 
681   Not Collective
682 
683   Input Parameters:
684 + n - number of values
685 . X - array of integers
686 . Y - second array of integers (first array of the pair)
687 - Z - third array of PetscCounts  (second array of the pair)
688 
689   Level: intermediate
690 
691   Note:
692   Usually X, Y are matrix row/column indices, and Z is a permutation array and therefore Z's type is PetscCount to allow 2B+ nonzeros even with 32-bit PetscInt.
693 
694 .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
695 @*/
696 PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
697 {
698   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
699   PetscCount t3;            /* temp for Z[] */
700 
701   PetscFunctionBegin;
702   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
703   PetscFunctionReturn(PETSC_SUCCESS);
704 }
705 
706 /*@
707   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
708 
709   Not Collective
710 
711   Input Parameters:
712 + n - number of values
713 - X - array of integers
714 
715   Output Parameter:
716 . sorted - flag whether the array is sorted
717 
718   Level: intermediate
719 
720 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
721 @*/
722 PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted)
723 {
724   PetscFunctionBegin;
725   PetscSorted(n, X, *sorted);
726   PetscFunctionReturn(PETSC_SUCCESS);
727 }
728 
729 /*@
730   PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
731 
732   Not Collective
733 
734   Input Parameters:
735 + n - number of values
736 - X - array of integers
737 
738   Level: intermediate
739 
740   Note:
741   This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
742   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
743   code to see which routine is fastest.
744 
745 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
746 @*/
747 PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[])
748 {
749   PetscMPIInt pivot, t1;
750 
751   PetscFunctionBegin;
752   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
753   PetscFunctionReturn(PETSC_SUCCESS);
754 }
755 
756 /*@
757   PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
758 
759   Not Collective
760 
761   Input Parameters:
762 + n - number of values
763 - X - array of integers
764 
765   Output Parameter:
766 . n - number of non-redundant values
767 
768   Level: intermediate
769 
770 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
771 @*/
772 PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
773 {
774   PetscInt s = 0, N = *n, b = 0;
775 
776   PetscFunctionBegin;
777   PetscCall(PetscSortMPIInt(N, X));
778   for (PetscInt i = 0; i < N - 1; i++) {
779     if (X[b + s + 1] != X[b]) {
780       X[b + 1] = X[b + s + 1];
781       b++;
782     } else s++;
783   }
784   *n = N - s;
785   PetscFunctionReturn(PETSC_SUCCESS);
786 }
787 
788 /*@
789   PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
790   changes a second `PetscMPIInt` array to match the sorted first array.
791 
792   Not Collective
793 
794   Input Parameters:
795 + n - number of values
796 . X - array of integers
797 - Y - second array of integers
798 
799   Level: intermediate
800 
801 .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
802 @*/
803 PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[])
804 {
805   PetscMPIInt pivot, t1, t2;
806 
807   PetscFunctionBegin;
808   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
809   PetscFunctionReturn(PETSC_SUCCESS);
810 }
811 
812 /*@
813   PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
814   changes a second array of `PetscInt` to match the sorted first array.
815 
816   Not Collective
817 
818   Input Parameters:
819 + n - number of values
820 . X - array of MPI integers
821 - Y - second array of Petsc integers
822 
823   Level: intermediate
824 
825   Note:
826   This routine is useful when one needs to sort MPI ranks with other integer arrays.
827 
828 .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
829 @*/
830 PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[])
831 {
832   PetscMPIInt pivot, t1;
833   PetscInt    t2;
834 
835   PetscFunctionBegin;
836   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
837   PetscFunctionReturn(PETSC_SUCCESS);
838 }
839 
840 /*@
841   PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
842   changes a second `PetscScalar` array to match the sorted first array.
843 
844   Not Collective
845 
846   Input Parameters:
847 + n - number of values
848 . X - array of integers
849 - Y - second array of scalars
850 
851   Level: intermediate
852 
853 .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
854 @*/
855 PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[])
856 {
857   PetscInt    pivot, t1;
858   PetscScalar t2;
859 
860   PetscFunctionBegin;
861   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
862   PetscFunctionReturn(PETSC_SUCCESS);
863 }
864 
865 /*@C
866   PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
867   changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
868   provide workspace (the size of an element in the data array) to use when sorting.
869 
870   Not Collective
871 
872   Input Parameters:
873 + n    - number of values
874 . X    - array of integers
875 . Y    - second array of data
876 . size - sizeof elements in the data array in bytes
877 - t2   - workspace of "size" bytes used when sorting
878 
879   Level: intermediate
880 
881 .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
882 @*/
883 PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2)
884 {
885   char    *YY = (char *)Y;
886   PetscInt t1, pivot, hi = n - 1;
887 
888   PetscFunctionBegin;
889   if (n < 8) {
890     for (PetscInt i = 0; i < n; i++) {
891       pivot = X[i];
892       for (PetscInt j = i + 1; j < n; j++) {
893         if (pivot > X[j]) {
894           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
895           pivot = X[i];
896         }
897       }
898     }
899   } else {
900     /* Two way partition */
901     PetscInt l = 0, r = hi;
902 
903     pivot = X[MEDIAN(X, hi)];
904     while (1) {
905       while (X[l] < pivot) l++;
906       while (X[r] > pivot) r--;
907       if (l >= r) {
908         r++;
909         break;
910       }
911       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
912       l++;
913       r--;
914     }
915     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
916     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
917   }
918   PetscFunctionReturn(PETSC_SUCCESS);
919 }
920 
921 /*@
922   PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
923 
924   Not Collective
925 
926   Input Parameters:
927 + an - number of values in the first array
928 . aI - first sorted array of integers
929 . bn - number of values in the second array
930 - bI - second array of integers
931 
932   Output Parameters:
933 + n - number of values in the merged array
934 - L - merged sorted array, this is allocated if an array is not provided
935 
936   Level: intermediate
937 
938 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
939 @*/
940 PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
941 {
942   PetscInt *L_ = *L, ak, bk, k;
943 
944   PetscFunctionBegin;
945   if (!L_) {
946     PetscCall(PetscMalloc1(an + bn, L));
947     L_ = *L;
948   }
949   k = ak = bk = 0;
950   while (ak < an && bk < bn) {
951     if (aI[ak] == bI[bk]) {
952       L_[k] = aI[ak];
953       ++ak;
954       ++bk;
955       ++k;
956     } else if (aI[ak] < bI[bk]) {
957       L_[k] = aI[ak];
958       ++ak;
959       ++k;
960     } else {
961       L_[k] = bI[bk];
962       ++bk;
963       ++k;
964     }
965   }
966   if (ak < an) {
967     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
968     k += (an - ak);
969   }
970   if (bk < bn) {
971     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
972     k += (bn - bk);
973   }
974   *n = k;
975   PetscFunctionReturn(PETSC_SUCCESS);
976 }
977 
978 /*@
979   PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
980   The additional arrays are the same length as sorted arrays and are merged
981   in the order determined by the merging of the sorted pair.
982 
983   Not Collective
984 
985   Input Parameters:
986 + an - number of values in the first array
987 . aI - first sorted array of integers
988 . aJ - first additional array of integers
989 . bn - number of values in the second array
990 . bI - second array of integers
991 - bJ - second additional of integers
992 
993   Output Parameters:
994 + n - number of values in the merged array (== an + bn)
995 . L - merged sorted array
996 - J - merged additional array
997 
998   Note:
999   if L or J point to non-null arrays then this routine will assume they are of the appropriate size and use them, otherwise this routine will allocate space for them
1000 
1001   Level: intermediate
1002 
1003 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1004 @*/
1005 PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J)
1006 {
1007   PetscInt n_, *L_, *J_, ak, bk, k;
1008 
1009   PetscFunctionBegin;
1010   PetscAssertPointer(L, 8);
1011   PetscAssertPointer(J, 9);
1012   n_ = an + bn;
1013   *n = n_;
1014   if (!*L) PetscCall(PetscMalloc1(n_, L));
1015   L_ = *L;
1016   if (!*J) PetscCall(PetscMalloc1(n_, J));
1017   J_ = *J;
1018   k = ak = bk = 0;
1019   while (ak < an && bk < bn) {
1020     if (aI[ak] <= bI[bk]) {
1021       L_[k] = aI[ak];
1022       J_[k] = aJ[ak];
1023       ++ak;
1024       ++k;
1025     } else {
1026       L_[k] = bI[bk];
1027       J_[k] = bJ[bk];
1028       ++bk;
1029       ++k;
1030     }
1031   }
1032   if (ak < an) {
1033     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
1034     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1035     k += (an - ak);
1036   }
1037   if (bk < bn) {
1038     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
1039     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1040   }
1041   PetscFunctionReturn(PETSC_SUCCESS);
1042 }
1043 
1044 /*@
1045   PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
1046 
1047   Not Collective
1048 
1049   Input Parameters:
1050 + an - number of values in the first array
1051 . aI - first sorted array of integers
1052 . bn - number of values in the second array
1053 - bI - second array of integers
1054 
1055   Output Parameters:
1056 + n - number of values in the merged array (<= an + bn)
1057 - L - merged sorted array, allocated if address of NULL pointer is passed
1058 
1059   Level: intermediate
1060 
1061 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1062 @*/
1063 PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1064 {
1065   PetscInt ai, bi, k;
1066 
1067   PetscFunctionBegin;
1068   if (!*L) PetscCall(PetscMalloc1((an + bn), L));
1069   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1070     PetscInt t = -1;
1071     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
1072     for (; bi < bn && bI[bi] == t; bi++)
1073       ;
1074     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
1075     for (; ai < an && aI[ai] == t; ai++)
1076       ;
1077   }
1078   *n = k;
1079   PetscFunctionReturn(PETSC_SUCCESS);
1080 }
1081 
1082 /*@C
1083   PetscProcessTree - Prepares tree data to be displayed graphically
1084 
1085   Not Collective
1086 
1087   Input Parameters:
1088 + n        - number of values
1089 . mask     - indicates those entries in the tree, location 0 is always masked
1090 - parentid - indicates the parent of each entry
1091 
1092   Output Parameters:
1093 + Nlevels   - the number of levels
1094 . Level     - for each node tells its level
1095 . Levelcnt  - the number of nodes on each level
1096 . Idbylevel - a list of ids on each of the levels, first level followed by second etc
1097 - Column    - for each id tells its column index
1098 
1099   Level: developer
1100 
1101   Note:
1102   This code is not currently used
1103 
1104 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
1105 @*/
1106 PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1107 {
1108   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1109   PetscBool done = PETSC_FALSE;
1110 
1111   PetscFunctionBegin;
1112   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
1113   for (i = 0; i < n; i++) {
1114     if (mask[i]) continue;
1115     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1116     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
1117   }
1118 
1119   for (i = 0; i < n; i++) {
1120     if (!mask[i]) nmask++;
1121   }
1122 
1123   /* determine the level in the tree of each node */
1124   PetscCall(PetscCalloc1(n, &level));
1125 
1126   level[0] = 1;
1127   while (!done) {
1128     done = PETSC_TRUE;
1129     for (i = 0; i < n; i++) {
1130       if (mask[i]) continue;
1131       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
1132       else if (!level[i]) done = PETSC_FALSE;
1133     }
1134   }
1135   for (i = 0; i < n; i++) {
1136     level[i]--;
1137     nlevels = PetscMax(nlevels, level[i]);
1138   }
1139 
1140   /* count the number of nodes on each level and its max */
1141   PetscCall(PetscCalloc1(nlevels, &levelcnt));
1142   for (i = 0; i < n; i++) {
1143     if (mask[i]) continue;
1144     levelcnt[level[i] - 1]++;
1145   }
1146   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
1147 
1148   /* for each level sort the ids by the parent id */
1149   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
1150   PetscCall(PetscMalloc1(nmask, &idbylevel));
1151   for (j = 1; j <= nlevels; j++) {
1152     cnt = 0;
1153     for (i = 0; i < n; i++) {
1154       if (mask[i]) continue;
1155       if (level[i] != j) continue;
1156       workid[cnt]         = i;
1157       workparentid[cnt++] = parentid[i];
1158     }
1159     /*  PetscIntView(cnt,workparentid,0);
1160     PetscIntView(cnt,workid,0);
1161     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
1162     PetscIntView(cnt,workparentid,0);
1163     PetscIntView(cnt,workid,0);*/
1164     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
1165     tcnt += cnt;
1166   }
1167   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
1168   PetscCall(PetscFree2(workid, workparentid));
1169 
1170   /* for each node list its column */
1171   PetscCall(PetscMalloc1(n, &column));
1172   cnt = 0;
1173   for (j = 0; j < nlevels; j++) {
1174     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
1175   }
1176 
1177   *Nlevels   = nlevels;
1178   *Level     = level;
1179   *Levelcnt  = levelcnt;
1180   *Idbylevel = idbylevel;
1181   *Column    = column;
1182   PetscFunctionReturn(PETSC_SUCCESS);
1183 }
1184 
1185 /*@
1186   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1187 
1188   Collective
1189 
1190   Input Parameters:
1191 + comm - the MPI communicator
1192 . n    - the local number of integers
1193 - keys - the local array of integers
1194 
1195   Output Parameters:
1196 . is_sorted - whether the array is globally sorted
1197 
1198   Level: developer
1199 
1200 .seealso: `PetscParallelSortInt()`
1201 @*/
1202 PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1203 {
1204   PetscBool   sorted;
1205   PetscInt    i, min, max, prevmax;
1206   PetscMPIInt rank;
1207 
1208   PetscFunctionBegin;
1209   sorted = PETSC_TRUE;
1210   min    = PETSC_MAX_INT;
1211   max    = PETSC_MIN_INT;
1212   if (n) {
1213     min = keys[0];
1214     max = keys[0];
1215   }
1216   for (i = 1; i < n; i++) {
1217     if (keys[i] < keys[i - 1]) break;
1218     min = PetscMin(min, keys[i]);
1219     max = PetscMax(max, keys[i]);
1220   }
1221   if (i < n) sorted = PETSC_FALSE;
1222   prevmax = PETSC_MIN_INT;
1223   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
1224   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1225   if (rank == 0) prevmax = PETSC_MIN_INT;
1226   if (prevmax > min) sorted = PETSC_FALSE;
1227   PetscCall(MPIU_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
1228   PetscFunctionReturn(PETSC_SUCCESS);
1229 }
1230