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