xref: /petsc/src/sys/utils/sorti.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
4c4762a1bSJed Brown    One can use src/sys/tests/ex52.c for benchmarking.
5e5c89e4eSSatish Balay  */
6af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"  I*/
7f1cab4e1SJunchao Zhang #include <petsc/private/hashseti.h>
8e5c89e4eSSatish Balay 
99371c9d4SSatish Balay #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))))
10a8582168SJed Brown 
11a8582168SJed Brown #define MEDIAN(v, right) MEDIAN3(v, right / 4, right / 2, right / 4 * 3)
12ef8e3583SJed Brown 
1359e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */
149371c9d4SSatish Balay #define SWAP1(a, b, t1) \
159371c9d4SSatish Balay   do { \
169371c9d4SSatish Balay     t1 = a; \
179371c9d4SSatish Balay     a  = b; \
189371c9d4SSatish Balay     b  = t1; \
199371c9d4SSatish Balay   } while (0)
209371c9d4SSatish Balay #define SWAP2(a, b, c, d, t1, t2) \
219371c9d4SSatish Balay   do { \
229371c9d4SSatish Balay     t1 = a; \
239371c9d4SSatish Balay     a  = b; \
249371c9d4SSatish Balay     b  = t1; \
259371c9d4SSatish Balay     t2 = c; \
269371c9d4SSatish Balay     c  = d; \
279371c9d4SSatish Balay     d  = t2; \
289371c9d4SSatish Balay   } while (0)
299371c9d4SSatish Balay #define SWAP3(a, b, c, d, e, f, t1, t2, t3) \
309371c9d4SSatish Balay   do { \
319371c9d4SSatish Balay     t1 = a; \
329371c9d4SSatish Balay     a  = b; \
339371c9d4SSatish Balay     b  = t1; \
349371c9d4SSatish Balay     t2 = c; \
359371c9d4SSatish Balay     c  = d; \
369371c9d4SSatish Balay     d  = t2; \
379371c9d4SSatish Balay     t3 = e; \
389371c9d4SSatish Balay     e  = f; \
399371c9d4SSatish Balay     f  = t3; \
409371c9d4SSatish Balay   } while (0)
4159e16d97SJunchao Zhang 
4259e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
4359e16d97SJunchao Zhang #define SWAP2Data(a, b, c, d, t1, t2, siz) \
4459e16d97SJunchao Zhang   do { \
459371c9d4SSatish Balay     t1 = a; \
469371c9d4SSatish Balay     a  = b; \
479371c9d4SSatish Balay     b  = t1; \
489566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(t2, c, siz)); \
499566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(c, d, siz)); \
509566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(d, t2, siz)); \
5159e16d97SJunchao Zhang   } while (0)
52e5c89e4eSSatish Balay 
53e5c89e4eSSatish Balay /*
5459e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
55e5c89e4eSSatish Balay 
5659e16d97SJunchao Zhang    Input Parameters:
5759e16d97SJunchao Zhang     + X         - array to partition
5859e16d97SJunchao Zhang     . pivot     - a pivot of X[]
5959e16d97SJunchao Zhang     . t1        - temp variable for X
6059e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
6159e16d97SJunchao Zhang 
6259e16d97SJunchao Zhang    Output Parameters:
6359e16d97SJunchao Zhang     + l,r       - of type PetscInt
6459e16d97SJunchao Zhang 
65811af0c4SBarry Smith    Note:
6659e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
6759e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
6859e16d97SJunchao Zhang  */
6959e16d97SJunchao Zhang #define TwoWayPartition1(X, pivot, t1, lo, hi, l, r) \
7059e16d97SJunchao Zhang   do { \
7159e16d97SJunchao Zhang     l = lo; \
7259e16d97SJunchao Zhang     r = hi; \
7359e16d97SJunchao Zhang     while (1) { \
7459e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
7559e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
769371c9d4SSatish Balay       if (l >= r) { \
779371c9d4SSatish Balay         r++; \
789371c9d4SSatish Balay         break; \
799371c9d4SSatish Balay       } \
8059e16d97SJunchao Zhang       SWAP1(X[l], X[r], t1); \
8159e16d97SJunchao Zhang       l++; \
8259e16d97SJunchao Zhang       r--; \
8359e16d97SJunchao Zhang     } \
8459e16d97SJunchao Zhang   } while (0)
8559e16d97SJunchao Zhang 
86ce605777SToby Isaac /*
87ce605777SToby Isaac    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
88ce605777SToby Isaac 
89ce605777SToby Isaac    Input Parameters:
90ce605777SToby Isaac     + X         - array to partition
91ce605777SToby Isaac     . pivot     - a pivot of X[]
92ce605777SToby Isaac     . t1        - temp variable for X
93ce605777SToby Isaac     - lo,hi     - lower and upper bound of the array
94ce605777SToby Isaac 
95ce605777SToby Isaac    Output Parameters:
96ce605777SToby Isaac     + l,r       - of type PetscInt
97ce605777SToby Isaac 
98811af0c4SBarry Smith    Note:
99ce605777SToby Isaac     The TwoWayPartition2/3 variants also partition other arrays along with X.
100ce605777SToby Isaac     These arrays can have different types, so they provide their own temp t2,t3
101ce605777SToby Isaac  */
102ce605777SToby Isaac #define TwoWayPartitionReverse1(X, pivot, t1, lo, hi, l, r) \
103ce605777SToby Isaac   do { \
104ce605777SToby Isaac     l = lo; \
105ce605777SToby Isaac     r = hi; \
106ce605777SToby Isaac     while (1) { \
107ce605777SToby Isaac       while (X[l] > pivot) l++; \
108ce605777SToby Isaac       while (X[r] < pivot) r--; \
1099371c9d4SSatish Balay       if (l >= r) { \
1109371c9d4SSatish Balay         r++; \
1119371c9d4SSatish Balay         break; \
1129371c9d4SSatish Balay       } \
113ce605777SToby Isaac       SWAP1(X[l], X[r], t1); \
114ce605777SToby Isaac       l++; \
115ce605777SToby Isaac       r--; \
116ce605777SToby Isaac     } \
117ce605777SToby Isaac   } while (0)
118ce605777SToby Isaac 
11959e16d97SJunchao Zhang #define TwoWayPartition2(X, Y, pivot, t1, t2, lo, hi, l, r) \
12059e16d97SJunchao Zhang   do { \
12159e16d97SJunchao Zhang     l = lo; \
12259e16d97SJunchao Zhang     r = hi; \
12359e16d97SJunchao Zhang     while (1) { \
12459e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
12559e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1269371c9d4SSatish Balay       if (l >= r) { \
1279371c9d4SSatish Balay         r++; \
1289371c9d4SSatish Balay         break; \
1299371c9d4SSatish Balay       } \
13059e16d97SJunchao Zhang       SWAP2(X[l], X[r], Y[l], Y[r], t1, t2); \
13159e16d97SJunchao Zhang       l++; \
13259e16d97SJunchao Zhang       r--; \
13359e16d97SJunchao Zhang     } \
13459e16d97SJunchao Zhang   } while (0)
13559e16d97SJunchao Zhang 
13659e16d97SJunchao Zhang #define TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, lo, hi, l, r) \
13759e16d97SJunchao Zhang   do { \
13859e16d97SJunchao Zhang     l = lo; \
13959e16d97SJunchao Zhang     r = hi; \
14059e16d97SJunchao Zhang     while (1) { \
14159e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
14259e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1439371c9d4SSatish Balay       if (l >= r) { \
1449371c9d4SSatish Balay         r++; \
1459371c9d4SSatish Balay         break; \
1469371c9d4SSatish Balay       } \
14759e16d97SJunchao Zhang       SWAP3(X[l], X[r], Y[l], Y[r], Z[l], Z[r], t1, t2, t3); \
14859e16d97SJunchao Zhang       l++; \
14959e16d97SJunchao Zhang       r--; \
15059e16d97SJunchao Zhang     } \
15159e16d97SJunchao Zhang   } while (0)
15259e16d97SJunchao Zhang 
15359e16d97SJunchao Zhang /* Templates for similar functions used below */
1545f80ce2aSJacob Faibussowitsch #define QuickSort1(FuncName, X, n, pivot, t1) \
15559e16d97SJunchao Zhang   do { \
156981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
15759e16d97SJunchao Zhang     if (n < 8) { \
15859e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
15959e16d97SJunchao Zhang         pivot = X[i]; \
16059e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
16159e16d97SJunchao Zhang           if (pivot > X[j]) { \
16259e16d97SJunchao Zhang             SWAP1(X[i], X[j], t1); \
16359e16d97SJunchao Zhang             pivot = X[i]; \
16459e16d97SJunchao Zhang           } \
16559e16d97SJunchao Zhang         } \
16659e16d97SJunchao Zhang       } \
16759e16d97SJunchao Zhang     } else { \
16859e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
16959e16d97SJunchao Zhang       pivot = X[p]; \
17059e16d97SJunchao Zhang       TwoWayPartition1(X, pivot, t1, 0, hi, l, r); \
1719566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1729566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
17359e16d97SJunchao Zhang     } \
17459e16d97SJunchao Zhang   } while (0)
17559e16d97SJunchao Zhang 
176ce605777SToby Isaac /* Templates for similar functions used below */
1775f80ce2aSJacob Faibussowitsch #define QuickSortReverse1(FuncName, X, n, pivot, t1) \
178ce605777SToby Isaac   do { \
179981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
180ce605777SToby Isaac     if (n < 8) { \
181ce605777SToby Isaac       for (i = 0; i < n; i++) { \
182ce605777SToby Isaac         pivot = X[i]; \
183ce605777SToby Isaac         for (j = i + 1; j < n; j++) { \
184ce605777SToby Isaac           if (pivot < X[j]) { \
185ce605777SToby Isaac             SWAP1(X[i], X[j], t1); \
186ce605777SToby Isaac             pivot = X[i]; \
187ce605777SToby Isaac           } \
188ce605777SToby Isaac         } \
189ce605777SToby Isaac       } \
190ce605777SToby Isaac     } else { \
191ce605777SToby Isaac       p     = MEDIAN(X, hi); \
192ce605777SToby Isaac       pivot = X[p]; \
193ce605777SToby Isaac       TwoWayPartitionReverse1(X, pivot, t1, 0, hi, l, r); \
1949566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1959566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
196ce605777SToby Isaac     } \
197ce605777SToby Isaac   } while (0)
198ce605777SToby Isaac 
1995f80ce2aSJacob Faibussowitsch #define QuickSort2(FuncName, X, Y, n, pivot, t1, t2) \
20059e16d97SJunchao Zhang   do { \
201981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
20259e16d97SJunchao Zhang     if (n < 8) { \
20359e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
20459e16d97SJunchao Zhang         pivot = X[i]; \
20559e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
20659e16d97SJunchao Zhang           if (pivot > X[j]) { \
20759e16d97SJunchao Zhang             SWAP2(X[i], X[j], Y[i], Y[j], t1, t2); \
20859e16d97SJunchao Zhang             pivot = X[i]; \
20959e16d97SJunchao Zhang           } \
21059e16d97SJunchao Zhang         } \
21159e16d97SJunchao Zhang       } \
21259e16d97SJunchao Zhang     } else { \
21359e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
21459e16d97SJunchao Zhang       pivot = X[p]; \
21559e16d97SJunchao Zhang       TwoWayPartition2(X, Y, pivot, t1, t2, 0, hi, l, r); \
2169566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y)); \
2179566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r)); \
21859e16d97SJunchao Zhang     } \
21959e16d97SJunchao Zhang   } while (0)
22059e16d97SJunchao Zhang 
2215f80ce2aSJacob Faibussowitsch #define QuickSort3(FuncName, X, Y, Z, n, pivot, t1, t2, t3) \
22259e16d97SJunchao Zhang   do { \
223981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
22459e16d97SJunchao Zhang     if (n < 8) { \
22559e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
22659e16d97SJunchao Zhang         pivot = X[i]; \
22759e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
22859e16d97SJunchao Zhang           if (pivot > X[j]) { \
22959e16d97SJunchao Zhang             SWAP3(X[i], X[j], Y[i], Y[j], Z[i], Z[j], t1, t2, t3); \
23059e16d97SJunchao Zhang             pivot = X[i]; \
23159e16d97SJunchao Zhang           } \
23259e16d97SJunchao Zhang         } \
23359e16d97SJunchao Zhang       } \
23459e16d97SJunchao Zhang     } else { \
23559e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
23659e16d97SJunchao Zhang       pivot = X[p]; \
23759e16d97SJunchao Zhang       TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, 0, hi, l, r); \
2389566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y, Z)); \
2399566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r, Z + r)); \
24059e16d97SJunchao Zhang     } \
24159e16d97SJunchao Zhang   } while (0)
242e5c89e4eSSatish Balay 
243e5c89e4eSSatish Balay /*@
244811af0c4SBarry Smith    PetscSortedInt - Determines whether the `PetscInt` array is sorted.
2456a7c706bSVaclav Hapla 
2466a7c706bSVaclav Hapla    Not Collective
2476a7c706bSVaclav Hapla 
2486a7c706bSVaclav Hapla    Input Parameters:
2496a7c706bSVaclav Hapla +  n  - number of values
2506a7c706bSVaclav Hapla -  X  - array of integers
2516a7c706bSVaclav Hapla 
252*2fe279fdSBarry Smith    Output Parameter:
2536a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
2546a7c706bSVaclav Hapla 
2556a7c706bSVaclav Hapla    Level: intermediate
2566a7c706bSVaclav Hapla 
257db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
2586a7c706bSVaclav Hapla @*/
259d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedInt(PetscInt n, const PetscInt X[], PetscBool *sorted)
260d71ae5a4SJacob Faibussowitsch {
2616a7c706bSVaclav Hapla   PetscFunctionBegin;
2625f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
2635f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(sorted, 3);
2646a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2666a7c706bSVaclav Hapla }
2676a7c706bSVaclav Hapla 
2686a7c706bSVaclav Hapla /*@
26962e5d2d2SJDBetteridge    PetscSortedInt64 - Determines whether the `PetscInt64` array is sorted.
27062e5d2d2SJDBetteridge 
27162e5d2d2SJDBetteridge    Not Collective
27262e5d2d2SJDBetteridge 
27362e5d2d2SJDBetteridge    Input Parameters:
27462e5d2d2SJDBetteridge +  n  - number of values
27562e5d2d2SJDBetteridge -  X  - array of integers
27662e5d2d2SJDBetteridge 
277*2fe279fdSBarry Smith    Output Parameter:
27862e5d2d2SJDBetteridge .  sorted - flag whether the array is sorted
27962e5d2d2SJDBetteridge 
28062e5d2d2SJDBetteridge    Level: intermediate
28162e5d2d2SJDBetteridge 
28262e5d2d2SJDBetteridge .seealso: `PetscSortInt64()`, `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
28362e5d2d2SJDBetteridge @*/
28462e5d2d2SJDBetteridge PetscErrorCode PetscSortedInt64(PetscInt n, const PetscInt64 X[], PetscBool *sorted)
28562e5d2d2SJDBetteridge {
28662e5d2d2SJDBetteridge   PetscFunctionBegin;
28762e5d2d2SJDBetteridge   if (n) PetscValidInt64Pointer(X, 2);
28862e5d2d2SJDBetteridge   PetscValidBoolPointer(sorted, 3);
28962e5d2d2SJDBetteridge   PetscSorted(n, X, *sorted);
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29162e5d2d2SJDBetteridge }
29262e5d2d2SJDBetteridge 
29362e5d2d2SJDBetteridge /*@
294811af0c4SBarry Smith    PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.
295e5c89e4eSSatish Balay 
296e5c89e4eSSatish Balay    Not Collective
297e5c89e4eSSatish Balay 
298e5c89e4eSSatish Balay    Input Parameters:
299e5c89e4eSSatish Balay +  n  - number of values
30059e16d97SJunchao Zhang -  X  - array of integers
301e5c89e4eSSatish Balay 
302811af0c4SBarry Smith    Note:
303811af0c4SBarry Smith    This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
304a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
305676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
306676f2a66SJacob Faibussowitsch 
307e5c89e4eSSatish Balay    Level: intermediate
308e5c89e4eSSatish Balay 
309db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
310e5c89e4eSSatish Balay @*/
311d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[])
312d71ae5a4SJacob Faibussowitsch {
31359e16d97SJunchao Zhang   PetscInt pivot, t1;
314e5c89e4eSSatish Balay 
315e5c89e4eSSatish Balay   PetscFunctionBegin;
3165f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
3175f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortInt, X, n, pivot, t1);
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319e5c89e4eSSatish Balay }
320e5c89e4eSSatish Balay 
3211db36a52SBarry Smith /*@
32262e5d2d2SJDBetteridge    PetscSortInt64 - Sorts an array of `PetscInt64` in place in increasing order.
32362e5d2d2SJDBetteridge 
32462e5d2d2SJDBetteridge    Not Collective
32562e5d2d2SJDBetteridge 
32662e5d2d2SJDBetteridge    Input Parameters:
32762e5d2d2SJDBetteridge +  n  - number of values
32862e5d2d2SJDBetteridge -  X  - array of integers
32962e5d2d2SJDBetteridge 
33062e5d2d2SJDBetteridge    Notes:
33162e5d2d2SJDBetteridge    This function sorts `PetscCount`s assumed to be in completely random order
33262e5d2d2SJDBetteridge 
33362e5d2d2SJDBetteridge    Level: intermediate
33462e5d2d2SJDBetteridge 
33562e5d2d2SJDBetteridge .seealso: `PetscSortInt()`
33662e5d2d2SJDBetteridge @*/
33762e5d2d2SJDBetteridge PetscErrorCode PetscSortInt64(PetscInt n, PetscInt64 X[])
33862e5d2d2SJDBetteridge {
33962e5d2d2SJDBetteridge   PetscCount pivot, t1;
34062e5d2d2SJDBetteridge 
34162e5d2d2SJDBetteridge   PetscFunctionBegin;
34262e5d2d2SJDBetteridge   if (n) PetscValidInt64Pointer(X, 2);
34362e5d2d2SJDBetteridge   QuickSort1(PetscSortInt64, X, n, pivot, t1);
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34562e5d2d2SJDBetteridge }
34662e5d2d2SJDBetteridge 
34762e5d2d2SJDBetteridge /*@
34862e5d2d2SJDBetteridge    PetscSortCount - Sorts an array of integers in place in increasing order.
34962e5d2d2SJDBetteridge 
35062e5d2d2SJDBetteridge    Not Collective
35162e5d2d2SJDBetteridge 
35262e5d2d2SJDBetteridge    Input Parameters:
35362e5d2d2SJDBetteridge +  n  - number of values
35462e5d2d2SJDBetteridge -  X  - array of integers
35562e5d2d2SJDBetteridge 
35662e5d2d2SJDBetteridge    Notes:
35762e5d2d2SJDBetteridge    This function sorts `PetscCount`s assumed to be in completely random order
35862e5d2d2SJDBetteridge 
35962e5d2d2SJDBetteridge    Level: intermediate
36062e5d2d2SJDBetteridge 
36162e5d2d2SJDBetteridge .seealso: `PetscSortInt()`
36262e5d2d2SJDBetteridge @*/
36362e5d2d2SJDBetteridge PetscErrorCode PetscSortCount(PetscInt n, PetscCount X[])
36462e5d2d2SJDBetteridge {
36562e5d2d2SJDBetteridge   PetscCount pivot, t1;
36662e5d2d2SJDBetteridge 
36762e5d2d2SJDBetteridge   PetscFunctionBegin;
36862e5d2d2SJDBetteridge   if (n) PetscValidCountPointer(X, 2);
36962e5d2d2SJDBetteridge   QuickSort1(PetscSortCount, X, n, pivot, t1);
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37162e5d2d2SJDBetteridge }
37262e5d2d2SJDBetteridge 
37362e5d2d2SJDBetteridge /*@
37462e5d2d2SJDBetteridge    PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
375ce605777SToby Isaac 
376ce605777SToby Isaac    Not Collective
377ce605777SToby Isaac 
378ce605777SToby Isaac    Input Parameters:
379ce605777SToby Isaac +  n  - number of values
380ce605777SToby Isaac -  X  - array of integers
381ce605777SToby Isaac 
382ce605777SToby Isaac    Level: intermediate
383ce605777SToby Isaac 
384db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
385ce605777SToby Isaac @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[])
387d71ae5a4SJacob Faibussowitsch {
388ce605777SToby Isaac   PetscInt pivot, t1;
389ce605777SToby Isaac 
390ce605777SToby Isaac   PetscFunctionBegin;
3915f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
3925f80ce2aSJacob Faibussowitsch   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
394ce605777SToby Isaac }
395ce605777SToby Isaac 
396ce605777SToby Isaac /*@
397811af0c4SBarry Smith    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array
39822ab5688SLisandro Dalcin 
39922ab5688SLisandro Dalcin    Not Collective
40022ab5688SLisandro Dalcin 
40122ab5688SLisandro Dalcin    Input Parameters:
40222ab5688SLisandro Dalcin +  n  - number of values
40359e16d97SJunchao Zhang -  X  - sorted array of integers
40422ab5688SLisandro Dalcin 
40522ab5688SLisandro Dalcin    Output Parameter:
40622ab5688SLisandro Dalcin .  n - number of non-redundant values
40722ab5688SLisandro Dalcin 
40822ab5688SLisandro Dalcin    Level: intermediate
40922ab5688SLisandro Dalcin 
410db781477SPatrick Sanan .seealso: `PetscSortInt()`
41122ab5688SLisandro Dalcin @*/
412d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[])
413d71ae5a4SJacob Faibussowitsch {
41422ab5688SLisandro Dalcin   PetscInt i, s = 0, N = *n, b = 0;
41522ab5688SLisandro Dalcin 
41622ab5688SLisandro Dalcin   PetscFunctionBegin;
4175f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n, 1);
418fb61b9e4SStefano Zampini   PetscCheckSorted(*n, X);
41922ab5688SLisandro Dalcin   for (i = 0; i < N - 1; i++) {
42059e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
4219371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
4229371c9d4SSatish Balay       b++;
42322ab5688SLisandro Dalcin     } else s++;
42422ab5688SLisandro Dalcin   }
42522ab5688SLisandro Dalcin   *n = N - s;
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42722ab5688SLisandro Dalcin }
42822ab5688SLisandro Dalcin 
42922ab5688SLisandro Dalcin /*@
430811af0c4SBarry Smith    PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates
431b6a18d21SVaclav Hapla 
432b6a18d21SVaclav Hapla    Not Collective
433b6a18d21SVaclav Hapla 
434b6a18d21SVaclav Hapla    Input Parameters:
435b6a18d21SVaclav Hapla +  n  - number of values
436b6a18d21SVaclav Hapla -  X  - sorted array of integers
437b6a18d21SVaclav Hapla 
438b6a18d21SVaclav Hapla    Output Parameter:
439b6a18d21SVaclav Hapla .  dups - True if the array has dups, otherwise false
440b6a18d21SVaclav Hapla 
441b6a18d21SVaclav Hapla    Level: intermediate
442b6a18d21SVaclav Hapla 
443db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
444b6a18d21SVaclav Hapla @*/
445d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg)
446d71ae5a4SJacob Faibussowitsch {
447b6a18d21SVaclav Hapla   PetscInt i;
448b6a18d21SVaclav Hapla 
449b6a18d21SVaclav Hapla   PetscFunctionBegin;
450b6a18d21SVaclav Hapla   PetscCheckSorted(n, X);
451b6a18d21SVaclav Hapla   *flg = PETSC_FALSE;
452b6a18d21SVaclav Hapla   for (i = 0; i < n - 1; i++) {
453b6a18d21SVaclav Hapla     if (X[i + 1] == X[i]) {
454b6a18d21SVaclav Hapla       *flg = PETSC_TRUE;
455b6a18d21SVaclav Hapla       break;
456b6a18d21SVaclav Hapla     }
457b6a18d21SVaclav Hapla   }
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459b6a18d21SVaclav Hapla }
460b6a18d21SVaclav Hapla 
461b6a18d21SVaclav Hapla /*@
462811af0c4SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
4631db36a52SBarry Smith 
4641db36a52SBarry Smith    Not Collective
4651db36a52SBarry Smith 
4661db36a52SBarry Smith    Input Parameters:
4671db36a52SBarry Smith +  n  - number of values
46859e16d97SJunchao Zhang -  X  - array of integers
4691db36a52SBarry Smith 
4701db36a52SBarry Smith    Output Parameter:
4711db36a52SBarry Smith .  n - number of non-redundant values
4721db36a52SBarry Smith 
4731db36a52SBarry Smith    Level: intermediate
4741db36a52SBarry Smith 
475db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
4761db36a52SBarry Smith @*/
477d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
478d71ae5a4SJacob Faibussowitsch {
4791db36a52SBarry Smith   PetscFunctionBegin;
4805f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n, 1);
4819566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(*n, X));
4829566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(n, X));
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4841db36a52SBarry Smith }
4851db36a52SBarry Smith 
48660e03357SMatthew G Knepley /*@
487811af0c4SBarry Smith  PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt`
48860e03357SMatthew G Knepley 
48960e03357SMatthew G Knepley    Not Collective
49060e03357SMatthew G Knepley 
49160e03357SMatthew G Knepley    Input Parameters:
49260e03357SMatthew G Knepley +  key - the integer to locate
493d9f1162dSJed Brown .  n   - number of values in the array
49459e16d97SJunchao Zhang -  X  - array of integers
49560e03357SMatthew G Knepley 
49660e03357SMatthew G Knepley    Output Parameter:
497d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
49860e03357SMatthew G Knepley 
49960e03357SMatthew G Knepley    Level: intermediate
50060e03357SMatthew G Knepley 
501db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
50260e03357SMatthew G Knepley @*/
503d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
504d71ae5a4SJacob Faibussowitsch {
505d9f1162dSJed Brown   PetscInt lo = 0, hi = n;
50660e03357SMatthew G Knepley 
50760e03357SMatthew G Knepley   PetscFunctionBegin;
5085f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(loc, 4);
5099371c9d4SSatish Balay   if (!n) {
5109371c9d4SSatish Balay     *loc = -1;
5113ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5129371c9d4SSatish Balay   }
513dadcf809SJacob Faibussowitsch   PetscValidIntPointer(X, 3);
5146a7c706bSVaclav Hapla   PetscCheckSorted(n, X);
515d9f1162dSJed Brown   while (hi - lo > 1) {
516d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo) / 2;
51759e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
518d9f1162dSJed Brown     else lo = mid;
51960e03357SMatthew G Knepley   }
52059e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52260e03357SMatthew G Knepley }
52360e03357SMatthew G Knepley 
524d2aeb606SJed Brown /*@
525811af0c4SBarry Smith   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
526f1cab4e1SJunchao Zhang 
527f1cab4e1SJunchao Zhang    Not Collective
528f1cab4e1SJunchao Zhang 
529f1cab4e1SJunchao Zhang    Input Parameters:
530f1cab4e1SJunchao Zhang +  n  - number of values in the array
531f1cab4e1SJunchao Zhang -  X  - array of integers
532f1cab4e1SJunchao Zhang 
533f1cab4e1SJunchao Zhang    Output Parameter:
534f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
535f1cab4e1SJunchao Zhang 
536f1cab4e1SJunchao Zhang    Level: intermediate
537f1cab4e1SJunchao Zhang 
538db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
539f1cab4e1SJunchao Zhang @*/
540d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
541d71ae5a4SJacob Faibussowitsch {
542f1cab4e1SJunchao Zhang   PetscInt   i;
543f1cab4e1SJunchao Zhang   PetscHSetI ht;
544f1cab4e1SJunchao Zhang   PetscBool  missing;
545f1cab4e1SJunchao Zhang 
546f1cab4e1SJunchao Zhang   PetscFunctionBegin;
5475f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
5485f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(dups, 3);
549f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
550f1cab4e1SJunchao Zhang   if (n > 1) {
5519566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&ht));
5529566063dSJacob Faibussowitsch     PetscCall(PetscHSetIResize(ht, n));
553f1cab4e1SJunchao Zhang     for (i = 0; i < n; i++) {
5549566063dSJacob Faibussowitsch       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
5559371c9d4SSatish Balay       if (!missing) {
5569371c9d4SSatish Balay         *dups = PETSC_TRUE;
5579371c9d4SSatish Balay         break;
5589371c9d4SSatish Balay       }
559f1cab4e1SJunchao Zhang     }
5609566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&ht));
561f1cab4e1SJunchao Zhang   }
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563f1cab4e1SJunchao Zhang }
564f1cab4e1SJunchao Zhang 
565f1cab4e1SJunchao Zhang /*@
566811af0c4SBarry Smith   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
567d2aeb606SJed Brown 
568d2aeb606SJed Brown    Not Collective
569d2aeb606SJed Brown 
570d2aeb606SJed Brown    Input Parameters:
571d2aeb606SJed Brown +  key - the integer to locate
572d2aeb606SJed Brown .  n   - number of values in the array
57359e16d97SJunchao Zhang -  X   - array of integers
574d2aeb606SJed Brown 
575d2aeb606SJed Brown    Output Parameter:
576d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
577d2aeb606SJed Brown 
578d2aeb606SJed Brown    Level: intermediate
579d2aeb606SJed Brown 
580db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
581d2aeb606SJed Brown @*/
582d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
583d71ae5a4SJacob Faibussowitsch {
584d2aeb606SJed Brown   PetscInt lo = 0, hi = n;
585d2aeb606SJed Brown 
586d2aeb606SJed Brown   PetscFunctionBegin;
587dadcf809SJacob Faibussowitsch   PetscValidIntPointer(loc, 4);
5889371c9d4SSatish Balay   if (!n) {
5899371c9d4SSatish Balay     *loc = -1;
5903ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5919371c9d4SSatish Balay   }
592dadcf809SJacob Faibussowitsch   PetscValidIntPointer(X, 3);
5936a7c706bSVaclav Hapla   PetscCheckSorted(n, X);
594d2aeb606SJed Brown   while (hi - lo > 1) {
595d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
59659e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
597d2aeb606SJed Brown     else lo = mid;
598d2aeb606SJed Brown   }
59959e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
6003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
601e5c89e4eSSatish Balay }
602e5c89e4eSSatish Balay 
603e5c89e4eSSatish Balay /*@
604811af0c4SBarry Smith    PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
605811af0c4SBarry Smith        changes a second array of `PetscInt` to match the sorted first array.
606e5c89e4eSSatish Balay 
607e5c89e4eSSatish Balay    Not Collective
608e5c89e4eSSatish Balay 
609e5c89e4eSSatish Balay    Input Parameters:
610e5c89e4eSSatish Balay +  n  - number of values
61159e16d97SJunchao Zhang .  X  - array of integers
61259e16d97SJunchao Zhang -  Y  - second array of integers
613e5c89e4eSSatish Balay 
614e5c89e4eSSatish Balay    Level: intermediate
615e5c89e4eSSatish Balay 
616db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
617e5c89e4eSSatish Balay @*/
618d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[])
619d71ae5a4SJacob Faibussowitsch {
62059e16d97SJunchao Zhang   PetscInt pivot, t1, t2;
621e5c89e4eSSatish Balay 
622e5c89e4eSSatish Balay   PetscFunctionBegin;
6235f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
6243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
625c1f0200aSDmitry Karpeev }
626c1f0200aSDmitry Karpeev 
627c1f0200aSDmitry Karpeev /*@
628811af0c4SBarry Smith    PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
629811af0c4SBarry Smith        changes a pair of `PetscInt` arrays to match the sorted first array.
630c1f0200aSDmitry Karpeev 
631c1f0200aSDmitry Karpeev    Not Collective
632c1f0200aSDmitry Karpeev 
633c1f0200aSDmitry Karpeev    Input Parameters:
634c1f0200aSDmitry Karpeev +  n  - number of values
63559e16d97SJunchao Zhang .  X  - array of integers
63659e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
63759e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
638c1f0200aSDmitry Karpeev 
639c1f0200aSDmitry Karpeev    Level: intermediate
640c1f0200aSDmitry Karpeev 
641db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
642c1f0200aSDmitry Karpeev @*/
643d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[])
644d71ae5a4SJacob Faibussowitsch {
64559e16d97SJunchao Zhang   PetscInt pivot, t1, t2, t3;
646c1f0200aSDmitry Karpeev 
647c1f0200aSDmitry Karpeev   PetscFunctionBegin;
6485f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650c1f0200aSDmitry Karpeev }
651c1f0200aSDmitry Karpeev 
65217d7d925SStefano Zampini /*@
653811af0c4SBarry Smith    PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
654811af0c4SBarry Smith        changes a second array of `PetscCount` to match the sorted first array.
655981bb840SJunchao Zhang 
656981bb840SJunchao Zhang    Not Collective
657981bb840SJunchao Zhang 
658981bb840SJunchao Zhang    Input Parameters:
659981bb840SJunchao Zhang +  n  - number of values
660981bb840SJunchao Zhang .  X  - array of integers
661981bb840SJunchao Zhang -  Y  - second array of PetscCounts (signed integers)
662981bb840SJunchao Zhang 
663981bb840SJunchao Zhang    Level: intermediate
664981bb840SJunchao Zhang 
665db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
666981bb840SJunchao Zhang @*/
667d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
668d71ae5a4SJacob Faibussowitsch {
669981bb840SJunchao Zhang   PetscInt   pivot, t1;
670981bb840SJunchao Zhang   PetscCount t2;
671981bb840SJunchao Zhang 
672981bb840SJunchao Zhang   PetscFunctionBegin;
6735f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
675981bb840SJunchao Zhang }
676981bb840SJunchao Zhang 
677981bb840SJunchao Zhang /*@
678811af0c4SBarry Smith    PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
679811af0c4SBarry Smith        changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
680981bb840SJunchao Zhang 
681981bb840SJunchao Zhang    Not Collective
682981bb840SJunchao Zhang 
683981bb840SJunchao Zhang    Input Parameters:
684981bb840SJunchao Zhang +  n  - number of values
685981bb840SJunchao Zhang .  X  - array of integers
686981bb840SJunchao Zhang .  Y  - second array of integers (first array of the pair)
687981bb840SJunchao Zhang -  Z  - third array of PetscCounts  (second array of the pair)
688981bb840SJunchao Zhang 
689981bb840SJunchao Zhang    Level: intermediate
690981bb840SJunchao Zhang 
691811af0c4SBarry Smith    Note:
692981bb840SJunchao Zhang     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.
693981bb840SJunchao Zhang 
694db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
695981bb840SJunchao Zhang @*/
696d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
697d71ae5a4SJacob Faibussowitsch {
698981bb840SJunchao Zhang   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
699981bb840SJunchao Zhang   PetscCount t3;            /* temp for Z[] */
700981bb840SJunchao Zhang 
701981bb840SJunchao Zhang   PetscFunctionBegin;
7025f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
704981bb840SJunchao Zhang }
705981bb840SJunchao Zhang 
706981bb840SJunchao Zhang /*@
707811af0c4SBarry Smith   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
7086a7c706bSVaclav Hapla 
7096a7c706bSVaclav Hapla    Not Collective
7106a7c706bSVaclav Hapla 
7116a7c706bSVaclav Hapla    Input Parameters:
7126a7c706bSVaclav Hapla +  n  - number of values
7136a7c706bSVaclav Hapla -  X  - array of integers
7146a7c706bSVaclav Hapla 
715*2fe279fdSBarry Smith    Output Parameter:
7166a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
7176a7c706bSVaclav Hapla 
7186a7c706bSVaclav Hapla    Level: intermediate
7196a7c706bSVaclav Hapla 
720db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
7216a7c706bSVaclav Hapla @*/
722d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted)
723d71ae5a4SJacob Faibussowitsch {
7246a7c706bSVaclav Hapla   PetscFunctionBegin;
7256a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7276a7c706bSVaclav Hapla }
7286a7c706bSVaclav Hapla 
7296a7c706bSVaclav Hapla /*@
730811af0c4SBarry Smith    PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
73117d7d925SStefano Zampini 
73217d7d925SStefano Zampini    Not Collective
73317d7d925SStefano Zampini 
73417d7d925SStefano Zampini    Input Parameters:
73517d7d925SStefano Zampini +  n  - number of values
73659e16d97SJunchao Zhang -  X  - array of integers
73717d7d925SStefano Zampini 
73817d7d925SStefano Zampini    Level: intermediate
73917d7d925SStefano Zampini 
740811af0c4SBarry Smith    Note:
741676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
742a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
743676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
744676f2a66SJacob Faibussowitsch 
745db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
74617d7d925SStefano Zampini @*/
747d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[])
748d71ae5a4SJacob Faibussowitsch {
74959e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
75017d7d925SStefano Zampini 
75117d7d925SStefano Zampini   PetscFunctionBegin;
7525f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
7533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75417d7d925SStefano Zampini }
75517d7d925SStefano Zampini 
75617d7d925SStefano Zampini /*@
757811af0c4SBarry Smith    PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
75817d7d925SStefano Zampini 
75917d7d925SStefano Zampini    Not Collective
76017d7d925SStefano Zampini 
76117d7d925SStefano Zampini    Input Parameters:
76217d7d925SStefano Zampini +  n  - number of values
76359e16d97SJunchao Zhang -  X  - array of integers
76417d7d925SStefano Zampini 
76517d7d925SStefano Zampini    Output Parameter:
76617d7d925SStefano Zampini .  n - number of non-redundant values
76717d7d925SStefano Zampini 
76817d7d925SStefano Zampini    Level: intermediate
76917d7d925SStefano Zampini 
770db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
77117d7d925SStefano Zampini @*/
772d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
773d71ae5a4SJacob Faibussowitsch {
7745f80ce2aSJacob Faibussowitsch   PetscInt s = 0, N = *n, b = 0;
77517d7d925SStefano Zampini 
77617d7d925SStefano Zampini   PetscFunctionBegin;
7779566063dSJacob Faibussowitsch   PetscCall(PetscSortMPIInt(N, X));
7785f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < N - 1; i++) {
77959e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
7809371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
7819371c9d4SSatish Balay       b++;
782a297a907SKarl Rupp     } else s++;
78317d7d925SStefano Zampini   }
78417d7d925SStefano Zampini   *n = N - s;
7853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78617d7d925SStefano Zampini }
787c1f0200aSDmitry Karpeev 
7884d615ea0SBarry Smith /*@
789811af0c4SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
790811af0c4SBarry Smith        changes a second `PetscMPIInt` array to match the sorted first array.
7914d615ea0SBarry Smith 
7924d615ea0SBarry Smith    Not Collective
7934d615ea0SBarry Smith 
7944d615ea0SBarry Smith    Input Parameters:
7954d615ea0SBarry Smith +  n  - number of values
79659e16d97SJunchao Zhang .  X  - array of integers
79759e16d97SJunchao Zhang -  Y  - second array of integers
7984d615ea0SBarry Smith 
7994d615ea0SBarry Smith    Level: intermediate
8004d615ea0SBarry Smith 
801db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
8024d615ea0SBarry Smith @*/
803d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[])
804d71ae5a4SJacob Faibussowitsch {
80559e16d97SJunchao Zhang   PetscMPIInt pivot, t1, t2;
8064d615ea0SBarry Smith 
8074d615ea0SBarry Smith   PetscFunctionBegin;
8085f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
8093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8105569a785SJunchao Zhang }
8115569a785SJunchao Zhang 
8125569a785SJunchao Zhang /*@
813811af0c4SBarry Smith    PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
814811af0c4SBarry Smith        changes a second array of `PetscInt` to match the sorted first array.
8155569a785SJunchao Zhang 
8165569a785SJunchao Zhang    Not Collective
8175569a785SJunchao Zhang 
8185569a785SJunchao Zhang    Input Parameters:
8195569a785SJunchao Zhang +  n  - number of values
82059e16d97SJunchao Zhang .  X  - array of MPI integers
82159e16d97SJunchao Zhang -  Y  - second array of Petsc integers
8225569a785SJunchao Zhang 
8235569a785SJunchao Zhang    Level: intermediate
8245569a785SJunchao Zhang 
825811af0c4SBarry Smith    Note:
826811af0c4SBarry Smith    This routine is useful when one needs to sort MPI ranks with other integer arrays.
8275569a785SJunchao Zhang 
828db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
8295569a785SJunchao Zhang @*/
830d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[])
831d71ae5a4SJacob Faibussowitsch {
83259e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
8335569a785SJunchao Zhang   PetscInt    t2;
8345569a785SJunchao Zhang 
8355569a785SJunchao Zhang   PetscFunctionBegin;
8365f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
8373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
838e5c89e4eSSatish Balay }
839e5c89e4eSSatish Balay 
840e5c89e4eSSatish Balay /*@
841811af0c4SBarry Smith    PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
842811af0c4SBarry Smith        changes a second `PetscScalar` array to match the sorted first array.
843e5c89e4eSSatish Balay 
844e5c89e4eSSatish Balay    Not Collective
845e5c89e4eSSatish Balay 
846e5c89e4eSSatish Balay    Input Parameters:
847e5c89e4eSSatish Balay +  n  - number of values
84859e16d97SJunchao Zhang .  X  - array of integers
84959e16d97SJunchao Zhang -  Y  - second array of scalars
850e5c89e4eSSatish Balay 
851e5c89e4eSSatish Balay    Level: intermediate
852e5c89e4eSSatish Balay 
853db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
854e5c89e4eSSatish Balay @*/
855d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[])
856d71ae5a4SJacob Faibussowitsch {
85759e16d97SJunchao Zhang   PetscInt    pivot, t1;
85859e16d97SJunchao Zhang   PetscScalar t2;
859e5c89e4eSSatish Balay 
860e5c89e4eSSatish Balay   PetscFunctionBegin;
8615f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
8623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86317df18f2SToby Isaac }
86417df18f2SToby Isaac 
8656c2863d0SBarry Smith /*@C
866811af0c4SBarry Smith    PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
86717df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
86817df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
86917df18f2SToby Isaac 
87017df18f2SToby Isaac    Not Collective
87117df18f2SToby Isaac 
87217df18f2SToby Isaac    Input Parameters:
87317df18f2SToby Isaac +  n  - number of values
87459e16d97SJunchao Zhang .  X  - array of integers
87559e16d97SJunchao Zhang .  Y  - second array of data
87617df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
87759e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
87817df18f2SToby Isaac 
87917df18f2SToby Isaac    Level: intermediate
88017df18f2SToby Isaac 
881db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
88217df18f2SToby Isaac @*/
883d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2)
884d71ae5a4SJacob Faibussowitsch {
88559e16d97SJunchao Zhang   char    *YY = (char *)Y;
8865f80ce2aSJacob Faibussowitsch   PetscInt t1, pivot, hi = n - 1;
88717df18f2SToby Isaac 
88817df18f2SToby Isaac   PetscFunctionBegin;
88917df18f2SToby Isaac   if (n < 8) {
8905f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < n; i++) {
89159e16d97SJunchao Zhang       pivot = X[i];
8925f80ce2aSJacob Faibussowitsch       for (PetscInt j = i + 1; j < n; j++) {
89359e16d97SJunchao Zhang         if (pivot > X[j]) {
89459e16d97SJunchao Zhang           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
89559e16d97SJunchao Zhang           pivot = X[i];
89617df18f2SToby Isaac         }
89717df18f2SToby Isaac       }
89817df18f2SToby Isaac     }
89917df18f2SToby Isaac   } else {
90059e16d97SJunchao Zhang     /* Two way partition */
9015f80ce2aSJacob Faibussowitsch     PetscInt l = 0, r = hi;
9025f80ce2aSJacob Faibussowitsch 
9035f80ce2aSJacob Faibussowitsch     pivot = X[MEDIAN(X, hi)];
90459e16d97SJunchao Zhang     while (1) {
90559e16d97SJunchao Zhang       while (X[l] < pivot) l++;
90659e16d97SJunchao Zhang       while (X[r] > pivot) r--;
9079371c9d4SSatish Balay       if (l >= r) {
9089371c9d4SSatish Balay         r++;
9099371c9d4SSatish Balay         break;
9109371c9d4SSatish Balay       }
91159e16d97SJunchao Zhang       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
91259e16d97SJunchao Zhang       l++;
91359e16d97SJunchao Zhang       r--;
91459e16d97SJunchao Zhang     }
9159566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
9169566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
91717df18f2SToby Isaac   }
9183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91917df18f2SToby Isaac }
92017df18f2SToby Isaac 
92121e72a00SBarry Smith /*@
922811af0c4SBarry Smith    PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
92321e72a00SBarry Smith 
92421e72a00SBarry Smith    Not Collective
92521e72a00SBarry Smith 
92621e72a00SBarry Smith    Input Parameters:
92721e72a00SBarry Smith +  an  - number of values in the first array
92821e72a00SBarry Smith .  aI  - first sorted array of integers
92921e72a00SBarry Smith .  bn  - number of values in the second array
93021e72a00SBarry Smith -  bI  - second array of integers
93121e72a00SBarry Smith 
93221e72a00SBarry Smith    Output Parameters:
93321e72a00SBarry Smith +  n   - number of values in the merged array
9346c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
93521e72a00SBarry Smith 
93621e72a00SBarry Smith    Level: intermediate
93721e72a00SBarry Smith 
938db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
93921e72a00SBarry Smith @*/
940d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
941d71ae5a4SJacob Faibussowitsch {
94221e72a00SBarry Smith   PetscInt *L_ = *L, ak, bk, k;
94321e72a00SBarry Smith 
944362febeeSStefano Zampini   PetscFunctionBegin;
94521e72a00SBarry Smith   if (!L_) {
9469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(an + bn, L));
94721e72a00SBarry Smith     L_ = *L;
94821e72a00SBarry Smith   }
94921e72a00SBarry Smith   k = ak = bk = 0;
95021e72a00SBarry Smith   while (ak < an && bk < bn) {
95121e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
95221e72a00SBarry Smith       L_[k] = aI[ak];
95321e72a00SBarry Smith       ++ak;
95421e72a00SBarry Smith       ++bk;
95521e72a00SBarry Smith       ++k;
95621e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
95721e72a00SBarry Smith       L_[k] = aI[ak];
95821e72a00SBarry Smith       ++ak;
95921e72a00SBarry Smith       ++k;
96021e72a00SBarry Smith     } else {
96121e72a00SBarry Smith       L_[k] = bI[bk];
96221e72a00SBarry Smith       ++bk;
96321e72a00SBarry Smith       ++k;
96421e72a00SBarry Smith     }
96521e72a00SBarry Smith   }
96621e72a00SBarry Smith   if (ak < an) {
9679566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
96821e72a00SBarry Smith     k += (an - ak);
96921e72a00SBarry Smith   }
97021e72a00SBarry Smith   if (bk < bn) {
9719566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
97221e72a00SBarry Smith     k += (bn - bk);
97321e72a00SBarry Smith   }
97421e72a00SBarry Smith   *n = k;
9753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97621e72a00SBarry Smith }
977b4301105SBarry Smith 
978c1f0200aSDmitry Karpeev /*@
979811af0c4SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
980c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
981c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
982c1f0200aSDmitry Karpeev 
983c1f0200aSDmitry Karpeev    Not Collective
984c1f0200aSDmitry Karpeev 
985c1f0200aSDmitry Karpeev    Input Parameters:
986c1f0200aSDmitry Karpeev +  an  - number of values in the first array
987c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
988c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
989c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
990c1f0200aSDmitry Karpeev .  bI  - second array of integers
991c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
992c1f0200aSDmitry Karpeev 
993c1f0200aSDmitry Karpeev    Output Parameters:
994c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
99514c49006SJed Brown .  L   - merged sorted array
996c1f0200aSDmitry Karpeev -  J   - merged additional array
997c1f0200aSDmitry Karpeev 
998811af0c4SBarry Smith    Note:
999da81f932SPierre Jolivet     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
1000811af0c4SBarry Smith 
1001c1f0200aSDmitry Karpeev    Level: intermediate
1002c1f0200aSDmitry Karpeev 
1003db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1004c1f0200aSDmitry Karpeev @*/
1005d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J)
1006d71ae5a4SJacob Faibussowitsch {
100770d8d27cSBarry Smith   PetscInt n_, *L_, *J_, ak, bk, k;
1008c1f0200aSDmitry Karpeev 
100914c49006SJed Brown   PetscFunctionBegin;
1010dadcf809SJacob Faibussowitsch   PetscValidPointer(L, 8);
1011dadcf809SJacob Faibussowitsch   PetscValidPointer(J, 9);
1012c1f0200aSDmitry Karpeev   n_ = an + bn;
1013c1f0200aSDmitry Karpeev   *n = n_;
101448a46eb9SPierre Jolivet   if (!*L) PetscCall(PetscMalloc1(n_, L));
1015d7aa01a8SHong Zhang   L_ = *L;
101648a46eb9SPierre Jolivet   if (!*J) PetscCall(PetscMalloc1(n_, J));
1017c1f0200aSDmitry Karpeev   J_ = *J;
1018c1f0200aSDmitry Karpeev   k = ak = bk = 0;
1019c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
1020c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
1021d7aa01a8SHong Zhang       L_[k] = aI[ak];
1022c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
1023c1f0200aSDmitry Karpeev       ++ak;
1024c1f0200aSDmitry Karpeev       ++k;
1025a297a907SKarl Rupp     } else {
1026d7aa01a8SHong Zhang       L_[k] = bI[bk];
1027c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
1028c1f0200aSDmitry Karpeev       ++bk;
1029c1f0200aSDmitry Karpeev       ++k;
1030c1f0200aSDmitry Karpeev     }
1031c1f0200aSDmitry Karpeev   }
1032c1f0200aSDmitry Karpeev   if (ak < an) {
10339566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
10349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1035c1f0200aSDmitry Karpeev     k += (an - ak);
1036c1f0200aSDmitry Karpeev   }
1037c1f0200aSDmitry Karpeev   if (bk < bn) {
10389566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
10399566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1040c1f0200aSDmitry Karpeev   }
10413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1042c1f0200aSDmitry Karpeev }
1043c1f0200aSDmitry Karpeev 
1044e498c390SJed Brown /*@
1045811af0c4SBarry Smith    PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
1046e498c390SJed Brown 
1047e498c390SJed Brown    Not Collective
1048e498c390SJed Brown 
1049e498c390SJed Brown    Input Parameters:
1050e498c390SJed Brown +  an  - number of values in the first array
1051e498c390SJed Brown .  aI  - first sorted array of integers
1052e498c390SJed Brown .  bn  - number of values in the second array
1053e498c390SJed Brown -  bI  - second array of integers
1054e498c390SJed Brown 
1055e498c390SJed Brown    Output Parameters:
1056e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
1057e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
1058e498c390SJed Brown 
1059e498c390SJed Brown    Level: intermediate
1060e498c390SJed Brown 
1061db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1062e498c390SJed Brown @*/
1063d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1064d71ae5a4SJacob Faibussowitsch {
1065e498c390SJed Brown   PetscInt ai, bi, k;
1066e498c390SJed Brown 
1067e498c390SJed Brown   PetscFunctionBegin;
10689566063dSJacob Faibussowitsch   if (!*L) PetscCall(PetscMalloc1((an + bn), L));
1069e498c390SJed Brown   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1070e498c390SJed Brown     PetscInt t = -1;
1071e498c390SJed Brown     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
10729371c9d4SSatish Balay     for (; bi < bn && bI[bi] == t; bi++)
10739371c9d4SSatish Balay       ;
1074e498c390SJed Brown     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
10759371c9d4SSatish Balay     for (; ai < an && aI[ai] == t; ai++)
10769371c9d4SSatish Balay       ;
1077e498c390SJed Brown   }
1078e498c390SJed Brown   *n = k;
10793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1080e498c390SJed Brown }
1081e5c89e4eSSatish Balay 
10826c2863d0SBarry Smith /*@C
108342eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
108442eaa7f3SBarry Smith 
108542eaa7f3SBarry Smith    Not Collective
108642eaa7f3SBarry Smith 
108742eaa7f3SBarry Smith    Input Parameters:
108842eaa7f3SBarry Smith +  n  - number of values
108942eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
109042eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
109142eaa7f3SBarry Smith 
109242eaa7f3SBarry Smith    Output Parameters:
109342eaa7f3SBarry Smith +  Nlevels - the number of levels
109442eaa7f3SBarry Smith .  Level - for each node tells its level
109542eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
109642eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
109742eaa7f3SBarry Smith -  Column - for each id tells its column index
109842eaa7f3SBarry Smith 
10996c2863d0SBarry Smith    Level: developer
110042eaa7f3SBarry Smith 
1101811af0c4SBarry Smith    Note:
110295452b02SPatrick Sanan     This code is not currently used
110342eaa7f3SBarry Smith 
1104db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
110542eaa7f3SBarry Smith @*/
1106d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1107d71ae5a4SJacob Faibussowitsch {
110842eaa7f3SBarry Smith   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1109ace3abfcSBarry Smith   PetscBool done = PETSC_FALSE;
111042eaa7f3SBarry Smith 
111142eaa7f3SBarry Smith   PetscFunctionBegin;
111208401ef6SPierre Jolivet   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
111342eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
111442eaa7f3SBarry Smith     if (mask[i]) continue;
111508401ef6SPierre Jolivet     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1116cc73adaaSBarry Smith     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
111742eaa7f3SBarry Smith   }
111842eaa7f3SBarry Smith 
111942eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
112042eaa7f3SBarry Smith     if (!mask[i]) nmask++;
112142eaa7f3SBarry Smith   }
112242eaa7f3SBarry Smith 
112342eaa7f3SBarry Smith   /* determine the level in the tree of each node */
11249566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &level));
1125a297a907SKarl Rupp 
112642eaa7f3SBarry Smith   level[0] = 1;
112742eaa7f3SBarry Smith   while (!done) {
112842eaa7f3SBarry Smith     done = PETSC_TRUE;
112942eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
113042eaa7f3SBarry Smith       if (mask[i]) continue;
113142eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
113242eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
113342eaa7f3SBarry Smith     }
113442eaa7f3SBarry Smith   }
113542eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
113642eaa7f3SBarry Smith     level[i]--;
113742eaa7f3SBarry Smith     nlevels = PetscMax(nlevels, level[i]);
113842eaa7f3SBarry Smith   }
113942eaa7f3SBarry Smith 
114042eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
11419566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nlevels, &levelcnt));
114242eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
114342eaa7f3SBarry Smith     if (mask[i]) continue;
114442eaa7f3SBarry Smith     levelcnt[level[i] - 1]++;
114542eaa7f3SBarry Smith   }
1146a297a907SKarl Rupp   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
114742eaa7f3SBarry Smith 
114842eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
11499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
11509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmask, &idbylevel));
115142eaa7f3SBarry Smith   for (j = 1; j <= nlevels; j++) {
115242eaa7f3SBarry Smith     cnt = 0;
115342eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
115442eaa7f3SBarry Smith       if (mask[i]) continue;
115542eaa7f3SBarry Smith       if (level[i] != j) continue;
115642eaa7f3SBarry Smith       workid[cnt]         = i;
115742eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
115842eaa7f3SBarry Smith     }
115942eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
116042eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
11619566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
116242eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
116342eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
11649566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
116542eaa7f3SBarry Smith     tcnt += cnt;
116642eaa7f3SBarry Smith   }
116708401ef6SPierre Jolivet   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
11689566063dSJacob Faibussowitsch   PetscCall(PetscFree2(workid, workparentid));
116942eaa7f3SBarry Smith 
117042eaa7f3SBarry Smith   /* for each node list its column */
11719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &column));
117242eaa7f3SBarry Smith   cnt = 0;
117342eaa7f3SBarry Smith   for (j = 0; j < nlevels; j++) {
1174ad540459SPierre Jolivet     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
117542eaa7f3SBarry Smith   }
117642eaa7f3SBarry Smith 
117742eaa7f3SBarry Smith   *Nlevels   = nlevels;
117842eaa7f3SBarry Smith   *Level     = level;
117942eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
118042eaa7f3SBarry Smith   *Idbylevel = idbylevel;
118142eaa7f3SBarry Smith   *Column    = column;
11823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118342eaa7f3SBarry Smith }
1184ce605777SToby Isaac 
1185ce605777SToby Isaac /*@
1186811af0c4SBarry Smith   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1187ce605777SToby Isaac 
1188ce605777SToby Isaac   Collective
1189ce605777SToby Isaac 
1190ce605777SToby Isaac   Input Parameters:
1191ce605777SToby Isaac + comm - the MPI communicator
1192ce605777SToby Isaac . n - the local number of integers
1193ce605777SToby Isaac - keys - the local array of integers
1194ce605777SToby Isaac 
1195ce605777SToby Isaac   Output Parameters:
1196ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1197ce605777SToby Isaac 
1198ce605777SToby Isaac   Level: developer
1199ce605777SToby Isaac 
1200db781477SPatrick Sanan .seealso: `PetscParallelSortInt()`
1201ce605777SToby Isaac @*/
1202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1203d71ae5a4SJacob Faibussowitsch {
1204ce605777SToby Isaac   PetscBool   sorted;
1205ce605777SToby Isaac   PetscInt    i, min, max, prevmax;
12062a1da528SToby Isaac   PetscMPIInt rank;
1207ce605777SToby Isaac 
1208ce605777SToby Isaac   PetscFunctionBegin;
1209ce605777SToby Isaac   sorted = PETSC_TRUE;
1210ce605777SToby Isaac   min    = PETSC_MAX_INT;
1211ce605777SToby Isaac   max    = PETSC_MIN_INT;
1212ce605777SToby Isaac   if (n) {
1213ce605777SToby Isaac     min = keys[0];
1214ce605777SToby Isaac     max = keys[0];
1215ce605777SToby Isaac   }
1216ce605777SToby Isaac   for (i = 1; i < n; i++) {
1217ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1218ce605777SToby Isaac     min = PetscMin(min, keys[i]);
1219ce605777SToby Isaac     max = PetscMax(max, keys[i]);
1220ce605777SToby Isaac   }
1221ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1222ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
12239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
12249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1225dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1226ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
12279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
12283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1229ce605777SToby Isaac }
1230