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