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