xref: /petsc/src/sys/utils/sorti.c (revision 811af0c4b09a35de4306c442f88bd09fdc09897d)
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 Parameters:
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   PetscFunctionBegin;
261   if (n) PetscValidIntPointer(X, 2);
262   PetscValidBoolPointer(sorted, 3);
263   PetscSorted(n, X, *sorted);
264   PetscFunctionReturn(0);
265 }
266 
267 /*@
268    PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.
269 
270    Not Collective
271 
272    Input Parameters:
273 +  n  - number of values
274 -  X  - array of integers
275 
276    Note:
277    This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
278    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
279    code to see which routine is fastest.
280 
281    Level: intermediate
282 
283 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
284 @*/
285 PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[]) {
286   PetscInt pivot, t1;
287 
288   PetscFunctionBegin;
289   if (n) PetscValidIntPointer(X, 2);
290   QuickSort1(PetscSortInt, X, n, pivot, t1);
291   PetscFunctionReturn(0);
292 }
293 
294 /*@
295    PetscSortReverseInt - Sorts an array of `PetscInt` in place in decreasing order.
296 
297    Not Collective
298 
299    Input Parameters:
300 +  n  - number of values
301 -  X  - array of integers
302 
303    Level: intermediate
304 
305 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
306 @*/
307 PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[]) {
308   PetscInt pivot, t1;
309 
310   PetscFunctionBegin;
311   if (n) PetscValidIntPointer(X, 2);
312   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
313   PetscFunctionReturn(0);
314 }
315 
316 /*@
317    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array
318 
319    Not Collective
320 
321    Input Parameters:
322 +  n  - number of values
323 -  X  - sorted array of integers
324 
325    Output Parameter:
326 .  n - number of non-redundant values
327 
328    Level: intermediate
329 
330 .seealso: `PetscSortInt()`
331 @*/
332 PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[]) {
333   PetscInt i, s = 0, N = *n, b = 0;
334 
335   PetscFunctionBegin;
336   PetscValidIntPointer(n, 1);
337   PetscCheckSorted(*n, X);
338   for (i = 0; i < N - 1; i++) {
339     if (X[b + s + 1] != X[b]) {
340       X[b + 1] = X[b + s + 1];
341       b++;
342     } else s++;
343   }
344   *n = N - s;
345   PetscFunctionReturn(0);
346 }
347 
348 /*@
349    PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates
350 
351    Not Collective
352 
353    Input Parameters:
354 +  n  - number of values
355 -  X  - sorted array of integers
356 
357    Output Parameter:
358 .  dups - True if the array has dups, otherwise false
359 
360    Level: intermediate
361 
362 .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
363 @*/
364 PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg) {
365   PetscInt i;
366 
367   PetscFunctionBegin;
368   PetscCheckSorted(n, X);
369   *flg = PETSC_FALSE;
370   for (i = 0; i < n - 1; i++) {
371     if (X[i + 1] == X[i]) {
372       *flg = PETSC_TRUE;
373       break;
374     }
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 /*@
380    PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
381 
382    Not Collective
383 
384    Input Parameters:
385 +  n  - number of values
386 -  X  - array of integers
387 
388    Output Parameter:
389 .  n - number of non-redundant values
390 
391    Level: intermediate
392 
393 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
394 @*/
395 PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[]) {
396   PetscFunctionBegin;
397   PetscValidIntPointer(n, 1);
398   PetscCall(PetscSortInt(*n, X));
399   PetscCall(PetscSortedRemoveDupsInt(n, X));
400   PetscFunctionReturn(0);
401 }
402 
403 /*@
404  PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt`
405 
406    Not Collective
407 
408    Input Parameters:
409 +  key - the integer to locate
410 .  n   - number of values in the array
411 -  X  - array of integers
412 
413    Output Parameter:
414 .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
415 
416    Level: intermediate
417 
418 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
419 @*/
420 PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc) {
421   PetscInt lo = 0, hi = n;
422 
423   PetscFunctionBegin;
424   PetscValidIntPointer(loc, 4);
425   if (!n) {
426     *loc = -1;
427     PetscFunctionReturn(0);
428   }
429   PetscValidIntPointer(X, 3);
430   PetscCheckSorted(n, X);
431   while (hi - lo > 1) {
432     PetscInt mid = lo + (hi - lo) / 2;
433     if (key < X[mid]) hi = mid;
434     else lo = mid;
435   }
436   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
437   PetscFunctionReturn(0);
438 }
439 
440 /*@
441   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
442 
443    Not Collective
444 
445    Input Parameters:
446 +  n  - number of values in the array
447 -  X  - array of integers
448 
449    Output Parameter:
450 .  dups - True if the array has dups, otherwise false
451 
452    Level: intermediate
453 
454 .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
455 @*/
456 PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups) {
457   PetscInt   i;
458   PetscHSetI ht;
459   PetscBool  missing;
460 
461   PetscFunctionBegin;
462   if (n) PetscValidIntPointer(X, 2);
463   PetscValidBoolPointer(dups, 3);
464   *dups = PETSC_FALSE;
465   if (n > 1) {
466     PetscCall(PetscHSetICreate(&ht));
467     PetscCall(PetscHSetIResize(ht, n));
468     for (i = 0; i < n; i++) {
469       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
470       if (!missing) {
471         *dups = PETSC_TRUE;
472         break;
473       }
474     }
475     PetscCall(PetscHSetIDestroy(&ht));
476   }
477   PetscFunctionReturn(0);
478 }
479 
480 /*@
481   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
482 
483    Not Collective
484 
485    Input Parameters:
486 +  key - the integer to locate
487 .  n   - number of values in the array
488 -  X   - array of integers
489 
490    Output Parameter:
491 .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
492 
493    Level: intermediate
494 
495 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
496 @*/
497 PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc) {
498   PetscInt lo = 0, hi = n;
499 
500   PetscFunctionBegin;
501   PetscValidIntPointer(loc, 4);
502   if (!n) {
503     *loc = -1;
504     PetscFunctionReturn(0);
505   }
506   PetscValidIntPointer(X, 3);
507   PetscCheckSorted(n, X);
508   while (hi - lo > 1) {
509     PetscInt mid = lo + (hi - lo) / 2;
510     if (key < X[mid]) hi = mid;
511     else lo = mid;
512   }
513   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
514   PetscFunctionReturn(0);
515 }
516 
517 /*@
518    PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
519        changes a second array of `PetscInt` to match the sorted first array.
520 
521    Not Collective
522 
523    Input Parameters:
524 +  n  - number of values
525 .  X  - array of integers
526 -  Y  - second array of integers
527 
528    Level: intermediate
529 
530 .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
531 @*/
532 PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[]) {
533   PetscInt pivot, t1, t2;
534 
535   PetscFunctionBegin;
536   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
537   PetscFunctionReturn(0);
538 }
539 
540 /*@
541    PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
542        changes a pair of `PetscInt` arrays to match the sorted first array.
543 
544    Not Collective
545 
546    Input Parameters:
547 +  n  - number of values
548 .  X  - array of integers
549 .  Y  - second array of integers (first array of the pair)
550 -  Z  - third array of integers  (second array of the pair)
551 
552    Level: intermediate
553 
554 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
555 @*/
556 PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[]) {
557   PetscInt pivot, t1, t2, t3;
558 
559   PetscFunctionBegin;
560   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
561   PetscFunctionReturn(0);
562 }
563 
564 /*@
565    PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
566        changes a second array of `PetscCount` to match the sorted first array.
567 
568    Not Collective
569 
570    Input Parameters:
571 +  n  - number of values
572 .  X  - array of integers
573 -  Y  - second array of PetscCounts (signed integers)
574 
575    Level: intermediate
576 
577 .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
578 @*/
579 PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[]) {
580   PetscInt   pivot, t1;
581   PetscCount t2;
582 
583   PetscFunctionBegin;
584   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
585   PetscFunctionReturn(0);
586 }
587 
588 /*@
589    PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
590        changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
591 
592    Not Collective
593 
594    Input Parameters:
595 +  n  - number of values
596 .  X  - array of integers
597 .  Y  - second array of integers (first array of the pair)
598 -  Z  - third array of PetscCounts  (second array of the pair)
599 
600    Level: intermediate
601 
602    Note:
603     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.
604 
605 .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
606 @*/
607 PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[]) {
608   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
609   PetscCount t3;            /* temp for Z[] */
610 
611   PetscFunctionBegin;
612   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
613   PetscFunctionReturn(0);
614 }
615 
616 /*@
617   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
618 
619    Not Collective
620 
621    Input Parameters:
622 +  n  - number of values
623 -  X  - array of integers
624 
625    Output Parameters:
626 .  sorted - flag whether the array is sorted
627 
628    Level: intermediate
629 
630 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
631 @*/
632 PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted) {
633   PetscFunctionBegin;
634   PetscSorted(n, X, *sorted);
635   PetscFunctionReturn(0);
636 }
637 
638 /*@
639    PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
640 
641    Not Collective
642 
643    Input Parameters:
644 +  n  - number of values
645 -  X  - array of integers
646 
647    Level: intermediate
648 
649    Note:
650    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
651    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
652    code to see which routine is fastest.
653 
654 .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
655 @*/
656 PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[]) {
657   PetscMPIInt pivot, t1;
658 
659   PetscFunctionBegin;
660   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
661   PetscFunctionReturn(0);
662 }
663 
664 /*@
665    PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
666 
667    Not Collective
668 
669    Input Parameters:
670 +  n  - number of values
671 -  X  - array of integers
672 
673    Output Parameter:
674 .  n - number of non-redundant values
675 
676    Level: intermediate
677 
678 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
679 @*/
680 PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[]) {
681   PetscInt s = 0, N = *n, b = 0;
682 
683   PetscFunctionBegin;
684   PetscCall(PetscSortMPIInt(N, X));
685   for (PetscInt i = 0; i < N - 1; i++) {
686     if (X[b + s + 1] != X[b]) {
687       X[b + 1] = X[b + s + 1];
688       b++;
689     } else s++;
690   }
691   *n = N - s;
692   PetscFunctionReturn(0);
693 }
694 
695 /*@
696    PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
697        changes a second `PetscMPIInt` array to match the sorted first array.
698 
699    Not Collective
700 
701    Input Parameters:
702 +  n  - number of values
703 .  X  - array of integers
704 -  Y  - second array of integers
705 
706    Level: intermediate
707 
708 .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
709 @*/
710 PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[]) {
711   PetscMPIInt pivot, t1, t2;
712 
713   PetscFunctionBegin;
714   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
715   PetscFunctionReturn(0);
716 }
717 
718 /*@
719    PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
720        changes a second array of `PetscInt` to match the sorted first array.
721 
722    Not Collective
723 
724    Input Parameters:
725 +  n  - number of values
726 .  X  - array of MPI integers
727 -  Y  - second array of Petsc integers
728 
729    Level: intermediate
730 
731    Note:
732    This routine is useful when one needs to sort MPI ranks with other integer arrays.
733 
734 .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
735 @*/
736 PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[]) {
737   PetscMPIInt pivot, t1;
738   PetscInt    t2;
739 
740   PetscFunctionBegin;
741   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
742   PetscFunctionReturn(0);
743 }
744 
745 /*@
746    PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
747        changes a second `PetscScalar` array to match the sorted first array.
748 
749    Not Collective
750 
751    Input Parameters:
752 +  n  - number of values
753 .  X  - array of integers
754 -  Y  - second array of scalars
755 
756    Level: intermediate
757 
758 .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
759 @*/
760 PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[]) {
761   PetscInt    pivot, t1;
762   PetscScalar t2;
763 
764   PetscFunctionBegin;
765   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
766   PetscFunctionReturn(0);
767 }
768 
769 /*@C
770    PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
771        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
772        provide workspace (the size of an element in the data array) to use when sorting.
773 
774    Not Collective
775 
776    Input Parameters:
777 +  n  - number of values
778 .  X  - array of integers
779 .  Y  - second array of data
780 .  size - sizeof elements in the data array in bytes
781 -  t2   - workspace of "size" bytes used when sorting
782 
783    Level: intermediate
784 
785 .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
786 @*/
787 PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2) {
788   char    *YY = (char *)Y;
789   PetscInt t1, pivot, hi = n - 1;
790 
791   PetscFunctionBegin;
792   if (n < 8) {
793     for (PetscInt i = 0; i < n; i++) {
794       pivot = X[i];
795       for (PetscInt j = i + 1; j < n; j++) {
796         if (pivot > X[j]) {
797           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
798           pivot = X[i];
799         }
800       }
801     }
802   } else {
803     /* Two way partition */
804     PetscInt l = 0, r = hi;
805 
806     pivot = X[MEDIAN(X, hi)];
807     while (1) {
808       while (X[l] < pivot) l++;
809       while (X[r] > pivot) r--;
810       if (l >= r) {
811         r++;
812         break;
813       }
814       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
815       l++;
816       r--;
817     }
818     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
819     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
820   }
821   PetscFunctionReturn(0);
822 }
823 
824 /*@
825    PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
826 
827    Not Collective
828 
829    Input Parameters:
830 +  an  - number of values in the first array
831 .  aI  - first sorted array of integers
832 .  bn  - number of values in the second array
833 -  bI  - second array of integers
834 
835    Output Parameters:
836 +  n   - number of values in the merged array
837 -  L   - merged sorted array, this is allocated if an array is not provided
838 
839    Level: intermediate
840 
841 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
842 @*/
843 PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) {
844   PetscInt *L_ = *L, ak, bk, k;
845 
846   PetscFunctionBegin;
847   if (!L_) {
848     PetscCall(PetscMalloc1(an + bn, L));
849     L_ = *L;
850   }
851   k = ak = bk = 0;
852   while (ak < an && bk < bn) {
853     if (aI[ak] == bI[bk]) {
854       L_[k] = aI[ak];
855       ++ak;
856       ++bk;
857       ++k;
858     } else if (aI[ak] < bI[bk]) {
859       L_[k] = aI[ak];
860       ++ak;
861       ++k;
862     } else {
863       L_[k] = bI[bk];
864       ++bk;
865       ++k;
866     }
867   }
868   if (ak < an) {
869     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
870     k += (an - ak);
871   }
872   if (bk < bn) {
873     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
874     k += (bn - bk);
875   }
876   *n = k;
877   PetscFunctionReturn(0);
878 }
879 
880 /*@
881    PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
882                                 The additional arrays are the same length as sorted arrays and are merged
883                                 in the order determined by the merging of the sorted pair.
884 
885    Not Collective
886 
887    Input Parameters:
888 +  an  - number of values in the first array
889 .  aI  - first sorted array of integers
890 .  aJ  - first additional array of integers
891 .  bn  - number of values in the second array
892 .  bI  - second array of integers
893 -  bJ  - second additional of integers
894 
895    Output Parameters:
896 +  n   - number of values in the merged array (== an + bn)
897 .  L   - merged sorted array
898 -  J   - merged additional array
899 
900    Note:
901     if L or J point to non-null arrays then this routine will assume they are of the approproate size and use them, otherwise this routine will allocate space for them
902 
903    Level: intermediate
904 
905 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
906 @*/
907 PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J) {
908   PetscInt n_, *L_, *J_, ak, bk, k;
909 
910   PetscFunctionBegin;
911   PetscValidPointer(L, 8);
912   PetscValidPointer(J, 9);
913   n_ = an + bn;
914   *n = n_;
915   if (!*L) PetscCall(PetscMalloc1(n_, L));
916   L_ = *L;
917   if (!*J) PetscCall(PetscMalloc1(n_, J));
918   J_ = *J;
919   k = ak = bk = 0;
920   while (ak < an && bk < bn) {
921     if (aI[ak] <= bI[bk]) {
922       L_[k] = aI[ak];
923       J_[k] = aJ[ak];
924       ++ak;
925       ++k;
926     } else {
927       L_[k] = bI[bk];
928       J_[k] = bJ[bk];
929       ++bk;
930       ++k;
931     }
932   }
933   if (ak < an) {
934     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
935     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
936     k += (an - ak);
937   }
938   if (bk < bn) {
939     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
940     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
941   }
942   PetscFunctionReturn(0);
943 }
944 
945 /*@
946    PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
947 
948    Not Collective
949 
950    Input Parameters:
951 +  an  - number of values in the first array
952 .  aI  - first sorted array of integers
953 .  bn  - number of values in the second array
954 -  bI  - second array of integers
955 
956    Output Parameters:
957 +  n   - number of values in the merged array (<= an + bn)
958 -  L   - merged sorted array, allocated if address of NULL pointer is passed
959 
960    Level: intermediate
961 
962 .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
963 @*/
964 PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L) {
965   PetscInt ai, bi, k;
966 
967   PetscFunctionBegin;
968   if (!*L) PetscCall(PetscMalloc1((an + bn), L));
969   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
970     PetscInt t = -1;
971     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
972     for (; bi < bn && bI[bi] == t; bi++)
973       ;
974     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
975     for (; ai < an && aI[ai] == t; ai++)
976       ;
977   }
978   *n = k;
979   PetscFunctionReturn(0);
980 }
981 
982 /*@C
983    PetscProcessTree - Prepares tree data to be displayed graphically
984 
985    Not Collective
986 
987    Input Parameters:
988 +  n  - number of values
989 .  mask - indicates those entries in the tree, location 0 is always masked
990 -  parentid - indicates the parent of each entry
991 
992    Output Parameters:
993 +  Nlevels - the number of levels
994 .  Level - for each node tells its level
995 .  Levelcnts - the number of nodes on each level
996 .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
997 -  Column - for each id tells its column index
998 
999    Level: developer
1000 
1001    Note:
1002     This code is not currently used
1003 
1004 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
1005 @*/
1006 PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column) {
1007   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1008   PetscBool done = PETSC_FALSE;
1009 
1010   PetscFunctionBegin;
1011   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
1012   for (i = 0; i < n; i++) {
1013     if (mask[i]) continue;
1014     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1015     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
1016   }
1017 
1018   for (i = 0; i < n; i++) {
1019     if (!mask[i]) nmask++;
1020   }
1021 
1022   /* determine the level in the tree of each node */
1023   PetscCall(PetscCalloc1(n, &level));
1024 
1025   level[0] = 1;
1026   while (!done) {
1027     done = PETSC_TRUE;
1028     for (i = 0; i < n; i++) {
1029       if (mask[i]) continue;
1030       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
1031       else if (!level[i]) done = PETSC_FALSE;
1032     }
1033   }
1034   for (i = 0; i < n; i++) {
1035     level[i]--;
1036     nlevels = PetscMax(nlevels, level[i]);
1037   }
1038 
1039   /* count the number of nodes on each level and its max */
1040   PetscCall(PetscCalloc1(nlevels, &levelcnt));
1041   for (i = 0; i < n; i++) {
1042     if (mask[i]) continue;
1043     levelcnt[level[i] - 1]++;
1044   }
1045   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
1046 
1047   /* for each level sort the ids by the parent id */
1048   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
1049   PetscCall(PetscMalloc1(nmask, &idbylevel));
1050   for (j = 1; j <= nlevels; j++) {
1051     cnt = 0;
1052     for (i = 0; i < n; i++) {
1053       if (mask[i]) continue;
1054       if (level[i] != j) continue;
1055       workid[cnt]         = i;
1056       workparentid[cnt++] = parentid[i];
1057     }
1058     /*  PetscIntView(cnt,workparentid,0);
1059     PetscIntView(cnt,workid,0);
1060     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
1061     PetscIntView(cnt,workparentid,0);
1062     PetscIntView(cnt,workid,0);*/
1063     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
1064     tcnt += cnt;
1065   }
1066   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
1067   PetscCall(PetscFree2(workid, workparentid));
1068 
1069   /* for each node list its column */
1070   PetscCall(PetscMalloc1(n, &column));
1071   cnt = 0;
1072   for (j = 0; j < nlevels; j++) {
1073     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
1074   }
1075 
1076   *Nlevels   = nlevels;
1077   *Level     = level;
1078   *Levelcnt  = levelcnt;
1079   *Idbylevel = idbylevel;
1080   *Column    = column;
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 /*@
1085   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1086 
1087   Collective
1088 
1089   Input Parameters:
1090 + comm - the MPI communicator
1091 . n - the local number of integers
1092 - keys - the local array of integers
1093 
1094   Output Parameters:
1095 . is_sorted - whether the array is globally sorted
1096 
1097   Level: developer
1098 
1099 .seealso: `PetscParallelSortInt()`
1100 @*/
1101 PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted) {
1102   PetscBool   sorted;
1103   PetscInt    i, min, max, prevmax;
1104   PetscMPIInt rank;
1105 
1106   PetscFunctionBegin;
1107   sorted = PETSC_TRUE;
1108   min    = PETSC_MAX_INT;
1109   max    = PETSC_MIN_INT;
1110   if (n) {
1111     min = keys[0];
1112     max = keys[0];
1113   }
1114   for (i = 1; i < n; i++) {
1115     if (keys[i] < keys[i - 1]) break;
1116     min = PetscMin(min, keys[i]);
1117     max = PetscMax(max, keys[i]);
1118   }
1119   if (i < n) sorted = PETSC_FALSE;
1120   prevmax = PETSC_MIN_INT;
1121   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
1122   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1123   if (rank == 0) prevmax = PETSC_MIN_INT;
1124   if (prevmax > min) sorted = PETSC_FALSE;
1125   PetscCallMPI(MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
1126   PetscFunctionReturn(0);
1127 }
1128