xref: /petsc/src/sys/utils/sorti.c (revision 42f965a0ff5eecf0e4aa5fc9705abe0e168bb37f)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
3c4762a1bSJed Brown    One can use src/sys/tests/ex52.c for benchmarking.
4e5c89e4eSSatish Balay  */
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"  I*/
6f1cab4e1SJunchao Zhang #include <petsc/private/hashseti.h>
7e5c89e4eSSatish Balay 
89371c9d4SSatish 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))))
9a8582168SJed Brown 
10a8582168SJed Brown #define MEDIAN(v, right) MEDIAN3(v, right / 4, right / 2, right / 4 * 3)
11ef8e3583SJed Brown 
1259e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */
139371c9d4SSatish Balay #define SWAP1(a, b, t1) \
149371c9d4SSatish Balay   do { \
159371c9d4SSatish Balay     t1 = a; \
169371c9d4SSatish Balay     a  = b; \
179371c9d4SSatish Balay     b  = t1; \
189371c9d4SSatish Balay   } while (0)
199371c9d4SSatish Balay #define SWAP2(a, b, c, d, t1, t2) \
209371c9d4SSatish Balay   do { \
219371c9d4SSatish Balay     t1 = a; \
229371c9d4SSatish Balay     a  = b; \
239371c9d4SSatish Balay     b  = t1; \
249371c9d4SSatish Balay     t2 = c; \
259371c9d4SSatish Balay     c  = d; \
269371c9d4SSatish Balay     d  = t2; \
279371c9d4SSatish Balay   } while (0)
289371c9d4SSatish Balay #define SWAP3(a, b, c, d, e, f, t1, t2, t3) \
299371c9d4SSatish Balay   do { \
309371c9d4SSatish Balay     t1 = a; \
319371c9d4SSatish Balay     a  = b; \
329371c9d4SSatish Balay     b  = t1; \
339371c9d4SSatish Balay     t2 = c; \
349371c9d4SSatish Balay     c  = d; \
359371c9d4SSatish Balay     d  = t2; \
369371c9d4SSatish Balay     t3 = e; \
379371c9d4SSatish Balay     e  = f; \
389371c9d4SSatish Balay     f  = t3; \
399371c9d4SSatish Balay   } while (0)
4059e16d97SJunchao Zhang 
4159e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
4259e16d97SJunchao Zhang #define SWAP2Data(a, b, c, d, t1, t2, siz) \
4359e16d97SJunchao Zhang   do { \
449371c9d4SSatish Balay     t1 = a; \
459371c9d4SSatish Balay     a  = b; \
469371c9d4SSatish Balay     b  = t1; \
479566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(t2, c, siz)); \
489566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(c, d, siz)); \
499566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(d, t2, siz)); \
5059e16d97SJunchao Zhang   } while (0)
51e5c89e4eSSatish Balay 
52e5c89e4eSSatish Balay /*
5359e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
54e5c89e4eSSatish Balay 
5559e16d97SJunchao Zhang    Input Parameters:
5659e16d97SJunchao Zhang     + X         - array to partition
5759e16d97SJunchao Zhang     . pivot     - a pivot of X[]
5859e16d97SJunchao Zhang     . t1        - temp variable for X
5959e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
6059e16d97SJunchao Zhang 
6159e16d97SJunchao Zhang    Output Parameters:
6259e16d97SJunchao Zhang     + l,r       - of type PetscInt
6359e16d97SJunchao Zhang 
64811af0c4SBarry Smith    Note:
6559e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
6659e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
6759e16d97SJunchao Zhang  */
6859e16d97SJunchao Zhang #define TwoWayPartition1(X, pivot, t1, lo, hi, l, r) \
6959e16d97SJunchao Zhang   do { \
7059e16d97SJunchao Zhang     l = lo; \
7159e16d97SJunchao Zhang     r = hi; \
7259e16d97SJunchao Zhang     while (1) { \
7359e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
7459e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
759371c9d4SSatish Balay       if (l >= r) { \
769371c9d4SSatish Balay         r++; \
779371c9d4SSatish Balay         break; \
789371c9d4SSatish Balay       } \
7959e16d97SJunchao Zhang       SWAP1(X[l], X[r], t1); \
8059e16d97SJunchao Zhang       l++; \
8159e16d97SJunchao Zhang       r--; \
8259e16d97SJunchao Zhang     } \
8359e16d97SJunchao Zhang   } while (0)
8459e16d97SJunchao Zhang 
85ce605777SToby Isaac /*
86ce605777SToby Isaac    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
87ce605777SToby Isaac 
88ce605777SToby Isaac    Input Parameters:
89ce605777SToby Isaac     + X         - array to partition
90ce605777SToby Isaac     . pivot     - a pivot of X[]
91ce605777SToby Isaac     . t1        - temp variable for X
92ce605777SToby Isaac     - lo,hi     - lower and upper bound of the array
93ce605777SToby Isaac 
94ce605777SToby Isaac    Output Parameters:
95ce605777SToby Isaac     + l,r       - of type PetscInt
96ce605777SToby Isaac 
97811af0c4SBarry Smith    Note:
98ce605777SToby Isaac     The TwoWayPartition2/3 variants also partition other arrays along with X.
99ce605777SToby Isaac     These arrays can have different types, so they provide their own temp t2,t3
100ce605777SToby Isaac  */
101ce605777SToby Isaac #define TwoWayPartitionReverse1(X, pivot, t1, lo, hi, l, r) \
102ce605777SToby Isaac   do { \
103ce605777SToby Isaac     l = lo; \
104ce605777SToby Isaac     r = hi; \
105ce605777SToby Isaac     while (1) { \
106ce605777SToby Isaac       while (X[l] > pivot) l++; \
107ce605777SToby Isaac       while (X[r] < pivot) r--; \
1089371c9d4SSatish Balay       if (l >= r) { \
1099371c9d4SSatish Balay         r++; \
1109371c9d4SSatish Balay         break; \
1119371c9d4SSatish Balay       } \
112ce605777SToby Isaac       SWAP1(X[l], X[r], t1); \
113ce605777SToby Isaac       l++; \
114ce605777SToby Isaac       r--; \
115ce605777SToby Isaac     } \
116ce605777SToby Isaac   } while (0)
117ce605777SToby Isaac 
11859e16d97SJunchao Zhang #define TwoWayPartition2(X, Y, pivot, t1, t2, lo, hi, l, r) \
11959e16d97SJunchao Zhang   do { \
12059e16d97SJunchao Zhang     l = lo; \
12159e16d97SJunchao Zhang     r = hi; \
12259e16d97SJunchao Zhang     while (1) { \
12359e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
12459e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1259371c9d4SSatish Balay       if (l >= r) { \
1269371c9d4SSatish Balay         r++; \
1279371c9d4SSatish Balay         break; \
1289371c9d4SSatish Balay       } \
12959e16d97SJunchao Zhang       SWAP2(X[l], X[r], Y[l], Y[r], t1, t2); \
13059e16d97SJunchao Zhang       l++; \
13159e16d97SJunchao Zhang       r--; \
13259e16d97SJunchao Zhang     } \
13359e16d97SJunchao Zhang   } while (0)
13459e16d97SJunchao Zhang 
13559e16d97SJunchao Zhang #define TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, lo, hi, l, r) \
13659e16d97SJunchao Zhang   do { \
13759e16d97SJunchao Zhang     l = lo; \
13859e16d97SJunchao Zhang     r = hi; \
13959e16d97SJunchao Zhang     while (1) { \
14059e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
14159e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1429371c9d4SSatish Balay       if (l >= r) { \
1439371c9d4SSatish Balay         r++; \
1449371c9d4SSatish Balay         break; \
1459371c9d4SSatish Balay       } \
14659e16d97SJunchao Zhang       SWAP3(X[l], X[r], Y[l], Y[r], Z[l], Z[r], t1, t2, t3); \
14759e16d97SJunchao Zhang       l++; \
14859e16d97SJunchao Zhang       r--; \
14959e16d97SJunchao Zhang     } \
15059e16d97SJunchao Zhang   } while (0)
15159e16d97SJunchao Zhang 
15259e16d97SJunchao Zhang /* Templates for similar functions used below */
1535f80ce2aSJacob Faibussowitsch #define QuickSort1(FuncName, X, n, pivot, t1) \
15459e16d97SJunchao Zhang   do { \
155981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
15659e16d97SJunchao Zhang     if (n < 8) { \
15759e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
15859e16d97SJunchao Zhang         pivot = X[i]; \
15959e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
16059e16d97SJunchao Zhang           if (pivot > X[j]) { \
16159e16d97SJunchao Zhang             SWAP1(X[i], X[j], t1); \
16259e16d97SJunchao Zhang             pivot = X[i]; \
16359e16d97SJunchao Zhang           } \
16459e16d97SJunchao Zhang         } \
16559e16d97SJunchao Zhang       } \
16659e16d97SJunchao Zhang     } else { \
16759e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
16859e16d97SJunchao Zhang       pivot = X[p]; \
16959e16d97SJunchao Zhang       TwoWayPartition1(X, pivot, t1, 0, hi, l, r); \
1709566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1719566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
17259e16d97SJunchao Zhang     } \
17359e16d97SJunchao Zhang   } while (0)
17459e16d97SJunchao Zhang 
175ce605777SToby Isaac /* Templates for similar functions used below */
1765f80ce2aSJacob Faibussowitsch #define QuickSortReverse1(FuncName, X, n, pivot, t1) \
177ce605777SToby Isaac   do { \
178981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
179ce605777SToby Isaac     if (n < 8) { \
180ce605777SToby Isaac       for (i = 0; i < n; i++) { \
181ce605777SToby Isaac         pivot = X[i]; \
182ce605777SToby Isaac         for (j = i + 1; j < n; j++) { \
183ce605777SToby Isaac           if (pivot < X[j]) { \
184ce605777SToby Isaac             SWAP1(X[i], X[j], t1); \
185ce605777SToby Isaac             pivot = X[i]; \
186ce605777SToby Isaac           } \
187ce605777SToby Isaac         } \
188ce605777SToby Isaac       } \
189ce605777SToby Isaac     } else { \
190ce605777SToby Isaac       p     = MEDIAN(X, hi); \
191ce605777SToby Isaac       pivot = X[p]; \
192ce605777SToby Isaac       TwoWayPartitionReverse1(X, pivot, t1, 0, hi, l, r); \
1939566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1949566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
195ce605777SToby Isaac     } \
196ce605777SToby Isaac   } while (0)
197ce605777SToby Isaac 
1985f80ce2aSJacob Faibussowitsch #define QuickSort2(FuncName, X, Y, n, pivot, t1, t2) \
19959e16d97SJunchao Zhang   do { \
200981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
20159e16d97SJunchao Zhang     if (n < 8) { \
20259e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
20359e16d97SJunchao Zhang         pivot = X[i]; \
20459e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
20559e16d97SJunchao Zhang           if (pivot > X[j]) { \
20659e16d97SJunchao Zhang             SWAP2(X[i], X[j], Y[i], Y[j], t1, t2); \
20759e16d97SJunchao Zhang             pivot = X[i]; \
20859e16d97SJunchao Zhang           } \
20959e16d97SJunchao Zhang         } \
21059e16d97SJunchao Zhang       } \
21159e16d97SJunchao Zhang     } else { \
21259e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
21359e16d97SJunchao Zhang       pivot = X[p]; \
21459e16d97SJunchao Zhang       TwoWayPartition2(X, Y, pivot, t1, t2, 0, hi, l, r); \
2159566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y)); \
2169566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r)); \
21759e16d97SJunchao Zhang     } \
21859e16d97SJunchao Zhang   } while (0)
21959e16d97SJunchao Zhang 
2205f80ce2aSJacob Faibussowitsch #define QuickSort3(FuncName, X, Y, Z, n, pivot, t1, t2, t3) \
22159e16d97SJunchao Zhang   do { \
222981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
22359e16d97SJunchao Zhang     if (n < 8) { \
22459e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
22559e16d97SJunchao Zhang         pivot = X[i]; \
22659e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
22759e16d97SJunchao Zhang           if (pivot > X[j]) { \
22859e16d97SJunchao Zhang             SWAP3(X[i], X[j], Y[i], Y[j], Z[i], Z[j], t1, t2, t3); \
22959e16d97SJunchao Zhang             pivot = X[i]; \
23059e16d97SJunchao Zhang           } \
23159e16d97SJunchao Zhang         } \
23259e16d97SJunchao Zhang       } \
23359e16d97SJunchao Zhang     } else { \
23459e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
23559e16d97SJunchao Zhang       pivot = X[p]; \
23659e16d97SJunchao Zhang       TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, 0, hi, l, r); \
2379566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y, Z)); \
2389566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r, Z + r)); \
23959e16d97SJunchao Zhang     } \
24059e16d97SJunchao Zhang   } while (0)
241e5c89e4eSSatish Balay 
242e5c89e4eSSatish Balay /*@
243811af0c4SBarry Smith   PetscSortedInt - Determines whether the `PetscInt` array is sorted.
2446a7c706bSVaclav Hapla 
2456a7c706bSVaclav Hapla   Not Collective
2466a7c706bSVaclav Hapla 
2476a7c706bSVaclav Hapla   Input Parameters:
2486a7c706bSVaclav Hapla + n - number of values
2496a7c706bSVaclav Hapla - X - array of integers
2506a7c706bSVaclav Hapla 
2512fe279fdSBarry Smith   Output Parameter:
2526a7c706bSVaclav Hapla . sorted - flag whether the array is sorted
2536a7c706bSVaclav Hapla 
2546a7c706bSVaclav Hapla   Level: intermediate
2556a7c706bSVaclav Hapla 
256db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
2576a7c706bSVaclav Hapla @*/
258d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedInt(PetscInt n, const PetscInt X[], PetscBool *sorted)
259d71ae5a4SJacob Faibussowitsch {
2606a7c706bSVaclav Hapla   PetscFunctionBegin;
2614f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
2624f572ea9SToby Isaac   PetscAssertPointer(sorted, 3);
2636a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2656a7c706bSVaclav Hapla }
2666a7c706bSVaclav Hapla 
2676a7c706bSVaclav Hapla /*@
26862e5d2d2SJDBetteridge   PetscSortedInt64 - Determines whether the `PetscInt64` array is sorted.
26962e5d2d2SJDBetteridge 
27062e5d2d2SJDBetteridge   Not Collective
27162e5d2d2SJDBetteridge 
27262e5d2d2SJDBetteridge   Input Parameters:
27362e5d2d2SJDBetteridge + n - number of values
27462e5d2d2SJDBetteridge - X - array of integers
27562e5d2d2SJDBetteridge 
2762fe279fdSBarry Smith   Output Parameter:
27762e5d2d2SJDBetteridge . sorted - flag whether the array is sorted
27862e5d2d2SJDBetteridge 
27962e5d2d2SJDBetteridge   Level: intermediate
28062e5d2d2SJDBetteridge 
28162e5d2d2SJDBetteridge .seealso: `PetscSortInt64()`, `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
28262e5d2d2SJDBetteridge @*/
28362e5d2d2SJDBetteridge PetscErrorCode PetscSortedInt64(PetscInt n, const PetscInt64 X[], PetscBool *sorted)
28462e5d2d2SJDBetteridge {
28562e5d2d2SJDBetteridge   PetscFunctionBegin;
2864f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
2874f572ea9SToby Isaac   PetscAssertPointer(sorted, 3);
28862e5d2d2SJDBetteridge   PetscSorted(n, X, *sorted);
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29062e5d2d2SJDBetteridge }
29162e5d2d2SJDBetteridge 
29262e5d2d2SJDBetteridge /*@
293811af0c4SBarry Smith   PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.
294e5c89e4eSSatish Balay 
295e5c89e4eSSatish Balay   Not Collective
296e5c89e4eSSatish Balay 
297e5c89e4eSSatish Balay   Input Parameters:
298e5c89e4eSSatish Balay + n - number of values
29959e16d97SJunchao Zhang - X - array of integers
300e5c89e4eSSatish Balay 
301811af0c4SBarry Smith   Note:
302811af0c4SBarry Smith   This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
303a5b23f4aSJose E. Roman   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
304676f2a66SJacob Faibussowitsch   code to see which routine is fastest.
305676f2a66SJacob Faibussowitsch 
306e5c89e4eSSatish Balay   Level: intermediate
307e5c89e4eSSatish Balay 
308db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
309e5c89e4eSSatish Balay @*/
310d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[])
311d71ae5a4SJacob Faibussowitsch {
31259e16d97SJunchao Zhang   PetscInt pivot, t1;
313e5c89e4eSSatish Balay 
314e5c89e4eSSatish Balay   PetscFunctionBegin;
3154f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
3165f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortInt, X, n, pivot, t1);
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318e5c89e4eSSatish Balay }
319e5c89e4eSSatish Balay 
3201db36a52SBarry Smith /*@
32162e5d2d2SJDBetteridge   PetscSortInt64 - Sorts an array of `PetscInt64` in place in increasing order.
32262e5d2d2SJDBetteridge 
32362e5d2d2SJDBetteridge   Not Collective
32462e5d2d2SJDBetteridge 
32562e5d2d2SJDBetteridge   Input Parameters:
32662e5d2d2SJDBetteridge + n - number of values
32762e5d2d2SJDBetteridge - X - array of integers
32862e5d2d2SJDBetteridge 
32962e5d2d2SJDBetteridge   Notes:
3307a43aa34SPierre Jolivet   This function sorts `PetscInt64`s assumed to be in completely random order
33162e5d2d2SJDBetteridge 
33262e5d2d2SJDBetteridge   Level: intermediate
33362e5d2d2SJDBetteridge 
33462e5d2d2SJDBetteridge .seealso: `PetscSortInt()`
33562e5d2d2SJDBetteridge @*/
33662e5d2d2SJDBetteridge PetscErrorCode PetscSortInt64(PetscInt n, PetscInt64 X[])
33762e5d2d2SJDBetteridge {
33862e5d2d2SJDBetteridge   PetscCount pivot, t1;
33962e5d2d2SJDBetteridge 
34062e5d2d2SJDBetteridge   PetscFunctionBegin;
3414f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
34262e5d2d2SJDBetteridge   QuickSort1(PetscSortInt64, X, n, pivot, t1);
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34462e5d2d2SJDBetteridge }
34562e5d2d2SJDBetteridge 
34662e5d2d2SJDBetteridge /*@
34762e5d2d2SJDBetteridge   PetscSortCount - Sorts an array of integers in place in increasing order.
34862e5d2d2SJDBetteridge 
34962e5d2d2SJDBetteridge   Not Collective
35062e5d2d2SJDBetteridge 
35162e5d2d2SJDBetteridge   Input Parameters:
35262e5d2d2SJDBetteridge + n - number of values
35362e5d2d2SJDBetteridge - X - array of integers
35462e5d2d2SJDBetteridge 
35562e5d2d2SJDBetteridge   Notes:
35662e5d2d2SJDBetteridge   This function sorts `PetscCount`s assumed to be in completely random order
35762e5d2d2SJDBetteridge 
35862e5d2d2SJDBetteridge   Level: intermediate
35962e5d2d2SJDBetteridge 
36062e5d2d2SJDBetteridge .seealso: `PetscSortInt()`
36162e5d2d2SJDBetteridge @*/
36262e5d2d2SJDBetteridge PetscErrorCode PetscSortCount(PetscInt n, PetscCount X[])
36362e5d2d2SJDBetteridge {
36462e5d2d2SJDBetteridge   PetscCount pivot, t1;
36562e5d2d2SJDBetteridge 
36662e5d2d2SJDBetteridge   PetscFunctionBegin;
3674f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
36862e5d2d2SJDBetteridge   QuickSort1(PetscSortCount, X, n, pivot, t1);
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37062e5d2d2SJDBetteridge }
37162e5d2d2SJDBetteridge 
37262e5d2d2SJDBetteridge /*@
37362e5d2d2SJDBetteridge   PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
374ce605777SToby Isaac 
375ce605777SToby Isaac   Not Collective
376ce605777SToby Isaac 
377ce605777SToby Isaac   Input Parameters:
378ce605777SToby Isaac + n - number of values
379ce605777SToby Isaac - X - array of integers
380ce605777SToby Isaac 
381ce605777SToby Isaac   Level: intermediate
382ce605777SToby Isaac 
383db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
384ce605777SToby Isaac @*/
385d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[])
386d71ae5a4SJacob Faibussowitsch {
387ce605777SToby Isaac   PetscInt pivot, t1;
388ce605777SToby Isaac 
389ce605777SToby Isaac   PetscFunctionBegin;
3904f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
3915f80ce2aSJacob Faibussowitsch   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393ce605777SToby Isaac }
394ce605777SToby Isaac 
395ce605777SToby Isaac /*@
396811af0c4SBarry Smith   PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array
39722ab5688SLisandro Dalcin 
39822ab5688SLisandro Dalcin   Not Collective
39922ab5688SLisandro Dalcin 
40022ab5688SLisandro Dalcin   Input Parameters:
40122ab5688SLisandro Dalcin + n - number of values
40259e16d97SJunchao Zhang - X - sorted array of integers
40322ab5688SLisandro Dalcin 
40422ab5688SLisandro Dalcin   Output Parameter:
40522ab5688SLisandro Dalcin . n - number of non-redundant values
40622ab5688SLisandro Dalcin 
40722ab5688SLisandro Dalcin   Level: intermediate
40822ab5688SLisandro Dalcin 
409db781477SPatrick Sanan .seealso: `PetscSortInt()`
41022ab5688SLisandro Dalcin @*/
411d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[])
412d71ae5a4SJacob Faibussowitsch {
41322ab5688SLisandro Dalcin   PetscInt i, s = 0, N = *n, b = 0;
41422ab5688SLisandro Dalcin 
41522ab5688SLisandro Dalcin   PetscFunctionBegin;
4164f572ea9SToby Isaac   PetscAssertPointer(n, 1);
417fb61b9e4SStefano Zampini   PetscCheckSorted(*n, X);
41822ab5688SLisandro Dalcin   for (i = 0; i < N - 1; i++) {
41959e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
4209371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
4219371c9d4SSatish Balay       b++;
42222ab5688SLisandro Dalcin     } else s++;
42322ab5688SLisandro Dalcin   }
42422ab5688SLisandro Dalcin   *n = N - s;
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42622ab5688SLisandro Dalcin }
42722ab5688SLisandro Dalcin 
42822ab5688SLisandro Dalcin /*@
429811af0c4SBarry Smith   PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates
430b6a18d21SVaclav Hapla 
431b6a18d21SVaclav Hapla   Not Collective
432b6a18d21SVaclav Hapla 
433b6a18d21SVaclav Hapla   Input Parameters:
434b6a18d21SVaclav Hapla + n - number of values
435b6a18d21SVaclav Hapla - X - sorted array of integers
436b6a18d21SVaclav Hapla 
437b6a18d21SVaclav Hapla   Output Parameter:
438aec76313SJacob Faibussowitsch . flg - True if the array has dups, otherwise false
439b6a18d21SVaclav Hapla 
440b6a18d21SVaclav Hapla   Level: intermediate
441b6a18d21SVaclav Hapla 
442db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
443b6a18d21SVaclav Hapla @*/
444d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg)
445d71ae5a4SJacob Faibussowitsch {
446b6a18d21SVaclav Hapla   PetscInt i;
447b6a18d21SVaclav Hapla 
448b6a18d21SVaclav Hapla   PetscFunctionBegin;
449b6a18d21SVaclav Hapla   PetscCheckSorted(n, X);
450b6a18d21SVaclav Hapla   *flg = PETSC_FALSE;
451b6a18d21SVaclav Hapla   for (i = 0; i < n - 1; i++) {
452b6a18d21SVaclav Hapla     if (X[i + 1] == X[i]) {
453b6a18d21SVaclav Hapla       *flg = PETSC_TRUE;
454b6a18d21SVaclav Hapla       break;
455b6a18d21SVaclav Hapla     }
456b6a18d21SVaclav Hapla   }
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458b6a18d21SVaclav Hapla }
459b6a18d21SVaclav Hapla 
460b6a18d21SVaclav Hapla /*@
461811af0c4SBarry Smith   PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
4621db36a52SBarry Smith 
4631db36a52SBarry Smith   Not Collective
4641db36a52SBarry Smith 
4651db36a52SBarry Smith   Input Parameters:
4661db36a52SBarry Smith + n - number of values
46759e16d97SJunchao Zhang - X - array of integers
4681db36a52SBarry Smith 
4691db36a52SBarry Smith   Output Parameter:
4701db36a52SBarry Smith . n - number of non-redundant values
4711db36a52SBarry Smith 
4721db36a52SBarry Smith   Level: intermediate
4731db36a52SBarry Smith 
474db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
4751db36a52SBarry Smith @*/
476d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
477d71ae5a4SJacob Faibussowitsch {
4781db36a52SBarry Smith   PetscFunctionBegin;
4794f572ea9SToby Isaac   PetscAssertPointer(n, 1);
4809566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(*n, X));
4819566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(n, X));
4823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4831db36a52SBarry Smith }
4841db36a52SBarry Smith 
48560e03357SMatthew G Knepley /*@
486811af0c4SBarry Smith   PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt`
48760e03357SMatthew G Knepley 
48860e03357SMatthew G Knepley   Not Collective
48960e03357SMatthew G Knepley 
49060e03357SMatthew G Knepley   Input Parameters:
49160e03357SMatthew G Knepley + key - the integer to locate
492d9f1162dSJed Brown . n   - number of values in the array
49359e16d97SJunchao Zhang - X   - array of integers
49460e03357SMatthew G Knepley 
49560e03357SMatthew G Knepley   Output Parameter:
496d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
49760e03357SMatthew G Knepley 
49860e03357SMatthew G Knepley   Level: intermediate
49960e03357SMatthew G Knepley 
500db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
50160e03357SMatthew G Knepley @*/
502d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
503d71ae5a4SJacob Faibussowitsch {
504d9f1162dSJed Brown   PetscInt lo = 0, hi = n;
50560e03357SMatthew G Knepley 
50660e03357SMatthew G Knepley   PetscFunctionBegin;
5074f572ea9SToby Isaac   PetscAssertPointer(loc, 4);
5089371c9d4SSatish Balay   if (!n) {
5099371c9d4SSatish Balay     *loc = -1;
5103ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5119371c9d4SSatish Balay   }
5124f572ea9SToby Isaac   PetscAssertPointer(X, 3);
513d9f1162dSJed Brown   while (hi - lo > 1) {
514d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo) / 2;
515*42f965a0SMatthew G. Knepley     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]);
51659e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
517d9f1162dSJed Brown     else lo = mid;
51860e03357SMatthew G Knepley   }
51959e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52160e03357SMatthew G Knepley }
52260e03357SMatthew G Knepley 
523d2aeb606SJed Brown /*@
524811af0c4SBarry Smith   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
525f1cab4e1SJunchao Zhang 
526f1cab4e1SJunchao Zhang   Not Collective
527f1cab4e1SJunchao Zhang 
528f1cab4e1SJunchao Zhang   Input Parameters:
529f1cab4e1SJunchao Zhang + n - number of values in the array
530f1cab4e1SJunchao Zhang - X - array of integers
531f1cab4e1SJunchao Zhang 
532f1cab4e1SJunchao Zhang   Output Parameter:
533f1cab4e1SJunchao Zhang . dups - True if the array has dups, otherwise false
534f1cab4e1SJunchao Zhang 
535f1cab4e1SJunchao Zhang   Level: intermediate
536f1cab4e1SJunchao Zhang 
537db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
538f1cab4e1SJunchao Zhang @*/
539d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
540d71ae5a4SJacob Faibussowitsch {
541f1cab4e1SJunchao Zhang   PetscInt   i;
542f1cab4e1SJunchao Zhang   PetscHSetI ht;
543f1cab4e1SJunchao Zhang   PetscBool  missing;
544f1cab4e1SJunchao Zhang 
545f1cab4e1SJunchao Zhang   PetscFunctionBegin;
5464f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
5474f572ea9SToby Isaac   PetscAssertPointer(dups, 3);
548f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
549f1cab4e1SJunchao Zhang   if (n > 1) {
5509566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&ht));
5519566063dSJacob Faibussowitsch     PetscCall(PetscHSetIResize(ht, n));
552f1cab4e1SJunchao Zhang     for (i = 0; i < n; i++) {
5539566063dSJacob Faibussowitsch       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
5549371c9d4SSatish Balay       if (!missing) {
5559371c9d4SSatish Balay         *dups = PETSC_TRUE;
5569371c9d4SSatish Balay         break;
5579371c9d4SSatish Balay       }
558f1cab4e1SJunchao Zhang     }
5599566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&ht));
560f1cab4e1SJunchao Zhang   }
5613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
562f1cab4e1SJunchao Zhang }
563f1cab4e1SJunchao Zhang 
564f1cab4e1SJunchao Zhang /*@
565811af0c4SBarry Smith   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
566d2aeb606SJed Brown 
567d2aeb606SJed Brown   Not Collective
568d2aeb606SJed Brown 
569d2aeb606SJed Brown   Input Parameters:
570d2aeb606SJed Brown + key - the integer to locate
571d2aeb606SJed Brown . n   - number of values in the array
57259e16d97SJunchao Zhang - X   - array of integers
573d2aeb606SJed Brown 
574d2aeb606SJed Brown   Output Parameter:
575d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
576d2aeb606SJed Brown 
577d2aeb606SJed Brown   Level: intermediate
578d2aeb606SJed Brown 
579db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
580d2aeb606SJed Brown @*/
581d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
582d71ae5a4SJacob Faibussowitsch {
583d2aeb606SJed Brown   PetscInt lo = 0, hi = n;
584d2aeb606SJed Brown 
585d2aeb606SJed Brown   PetscFunctionBegin;
5864f572ea9SToby Isaac   PetscAssertPointer(loc, 4);
5879371c9d4SSatish Balay   if (!n) {
5889371c9d4SSatish Balay     *loc = -1;
5893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5909371c9d4SSatish Balay   }
5914f572ea9SToby Isaac   PetscAssertPointer(X, 3);
592d2aeb606SJed Brown   while (hi - lo > 1) {
593d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
594*42f965a0SMatthew G. Knepley     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]);
59559e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
596d2aeb606SJed Brown     else lo = mid;
597d2aeb606SJed Brown   }
59859e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
5993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
600e5c89e4eSSatish Balay }
601e5c89e4eSSatish Balay 
602e5c89e4eSSatish Balay /*@
603811af0c4SBarry Smith   PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
604811af0c4SBarry Smith   changes a second array of `PetscInt` to match the sorted first array.
605e5c89e4eSSatish Balay 
606e5c89e4eSSatish Balay   Not Collective
607e5c89e4eSSatish Balay 
608e5c89e4eSSatish Balay   Input Parameters:
609e5c89e4eSSatish Balay + n - number of values
61059e16d97SJunchao Zhang . X - array of integers
61159e16d97SJunchao Zhang - Y - second array of integers
612e5c89e4eSSatish Balay 
613e5c89e4eSSatish Balay   Level: intermediate
614e5c89e4eSSatish Balay 
615db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
616e5c89e4eSSatish Balay @*/
617d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[])
618d71ae5a4SJacob Faibussowitsch {
61959e16d97SJunchao Zhang   PetscInt pivot, t1, t2;
620e5c89e4eSSatish Balay 
621e5c89e4eSSatish Balay   PetscFunctionBegin;
6225f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
6233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
624c1f0200aSDmitry Karpeev }
625c1f0200aSDmitry Karpeev 
626c1f0200aSDmitry Karpeev /*@
627811af0c4SBarry Smith   PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
628811af0c4SBarry Smith   changes a pair of `PetscInt` arrays to match the sorted first array.
629c1f0200aSDmitry Karpeev 
630c1f0200aSDmitry Karpeev   Not Collective
631c1f0200aSDmitry Karpeev 
632c1f0200aSDmitry Karpeev   Input Parameters:
633c1f0200aSDmitry Karpeev + n - number of values
63459e16d97SJunchao Zhang . X - array of integers
63559e16d97SJunchao Zhang . Y - second array of integers (first array of the pair)
63659e16d97SJunchao Zhang - Z - third array of integers  (second array of the pair)
637c1f0200aSDmitry Karpeev 
638c1f0200aSDmitry Karpeev   Level: intermediate
639c1f0200aSDmitry Karpeev 
640db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
641c1f0200aSDmitry Karpeev @*/
642d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[])
643d71ae5a4SJacob Faibussowitsch {
64459e16d97SJunchao Zhang   PetscInt pivot, t1, t2, t3;
645c1f0200aSDmitry Karpeev 
646c1f0200aSDmitry Karpeev   PetscFunctionBegin;
6475f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
649c1f0200aSDmitry Karpeev }
650c1f0200aSDmitry Karpeev 
65117d7d925SStefano Zampini /*@
652811af0c4SBarry Smith   PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
653811af0c4SBarry Smith   changes a second array of `PetscCount` to match the sorted first array.
654981bb840SJunchao Zhang 
655981bb840SJunchao Zhang   Not Collective
656981bb840SJunchao Zhang 
657981bb840SJunchao Zhang   Input Parameters:
658981bb840SJunchao Zhang + n - number of values
659981bb840SJunchao Zhang . X - array of integers
660981bb840SJunchao Zhang - Y - second array of PetscCounts (signed integers)
661981bb840SJunchao Zhang 
662981bb840SJunchao Zhang   Level: intermediate
663981bb840SJunchao Zhang 
664db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
665981bb840SJunchao Zhang @*/
666d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
667d71ae5a4SJacob Faibussowitsch {
668981bb840SJunchao Zhang   PetscInt   pivot, t1;
669981bb840SJunchao Zhang   PetscCount t2;
670981bb840SJunchao Zhang 
671981bb840SJunchao Zhang   PetscFunctionBegin;
6725f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
6733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
674981bb840SJunchao Zhang }
675981bb840SJunchao Zhang 
676981bb840SJunchao Zhang /*@
677811af0c4SBarry Smith   PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
678811af0c4SBarry Smith   changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
679981bb840SJunchao Zhang 
680981bb840SJunchao Zhang   Not Collective
681981bb840SJunchao Zhang 
682981bb840SJunchao Zhang   Input Parameters:
683981bb840SJunchao Zhang + n - number of values
684981bb840SJunchao Zhang . X - array of integers
685981bb840SJunchao Zhang . Y - second array of integers (first array of the pair)
686981bb840SJunchao Zhang - Z - third array of PetscCounts  (second array of the pair)
687981bb840SJunchao Zhang 
688981bb840SJunchao Zhang   Level: intermediate
689981bb840SJunchao Zhang 
690811af0c4SBarry Smith   Note:
691981bb840SJunchao 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.
692981bb840SJunchao Zhang 
693db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
694981bb840SJunchao Zhang @*/
695d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
696d71ae5a4SJacob Faibussowitsch {
697981bb840SJunchao Zhang   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
698981bb840SJunchao Zhang   PetscCount t3;            /* temp for Z[] */
699981bb840SJunchao Zhang 
700981bb840SJunchao Zhang   PetscFunctionBegin;
7015f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
703981bb840SJunchao Zhang }
704981bb840SJunchao Zhang 
705981bb840SJunchao Zhang /*@
706811af0c4SBarry Smith   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
7076a7c706bSVaclav Hapla 
7086a7c706bSVaclav Hapla   Not Collective
7096a7c706bSVaclav Hapla 
7106a7c706bSVaclav Hapla   Input Parameters:
7116a7c706bSVaclav Hapla + n - number of values
7126a7c706bSVaclav Hapla - X - array of integers
7136a7c706bSVaclav Hapla 
7142fe279fdSBarry Smith   Output Parameter:
7156a7c706bSVaclav Hapla . sorted - flag whether the array is sorted
7166a7c706bSVaclav Hapla 
7176a7c706bSVaclav Hapla   Level: intermediate
7186a7c706bSVaclav Hapla 
719db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
7206a7c706bSVaclav Hapla @*/
721d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted)
722d71ae5a4SJacob Faibussowitsch {
7236a7c706bSVaclav Hapla   PetscFunctionBegin;
7246a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
7253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7266a7c706bSVaclav Hapla }
7276a7c706bSVaclav Hapla 
7286a7c706bSVaclav Hapla /*@
729811af0c4SBarry Smith   PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
73017d7d925SStefano Zampini 
73117d7d925SStefano Zampini   Not Collective
73217d7d925SStefano Zampini 
73317d7d925SStefano Zampini   Input Parameters:
73417d7d925SStefano Zampini + n - number of values
73559e16d97SJunchao Zhang - X - array of integers
73617d7d925SStefano Zampini 
73717d7d925SStefano Zampini   Level: intermediate
73817d7d925SStefano Zampini 
739811af0c4SBarry Smith   Note:
740676f2a66SJacob Faibussowitsch   This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
741a5b23f4aSJose E. Roman   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
742676f2a66SJacob Faibussowitsch   code to see which routine is fastest.
743676f2a66SJacob Faibussowitsch 
744db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
74517d7d925SStefano Zampini @*/
746d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[])
747d71ae5a4SJacob Faibussowitsch {
74859e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
74917d7d925SStefano Zampini 
75017d7d925SStefano Zampini   PetscFunctionBegin;
7515f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
7523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75317d7d925SStefano Zampini }
75417d7d925SStefano Zampini 
75517d7d925SStefano Zampini /*@
756811af0c4SBarry Smith   PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
75717d7d925SStefano Zampini 
75817d7d925SStefano Zampini   Not Collective
75917d7d925SStefano Zampini 
76017d7d925SStefano Zampini   Input Parameters:
76117d7d925SStefano Zampini + n - number of values
76259e16d97SJunchao Zhang - X - array of integers
76317d7d925SStefano Zampini 
76417d7d925SStefano Zampini   Output Parameter:
76517d7d925SStefano Zampini . n - number of non-redundant values
76617d7d925SStefano Zampini 
76717d7d925SStefano Zampini   Level: intermediate
76817d7d925SStefano Zampini 
769db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
77017d7d925SStefano Zampini @*/
771d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
772d71ae5a4SJacob Faibussowitsch {
7735f80ce2aSJacob Faibussowitsch   PetscInt s = 0, N = *n, b = 0;
77417d7d925SStefano Zampini 
77517d7d925SStefano Zampini   PetscFunctionBegin;
7769566063dSJacob Faibussowitsch   PetscCall(PetscSortMPIInt(N, X));
7775f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < N - 1; i++) {
77859e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
7799371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
7809371c9d4SSatish Balay       b++;
781a297a907SKarl Rupp     } else s++;
78217d7d925SStefano Zampini   }
78317d7d925SStefano Zampini   *n = N - s;
7843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78517d7d925SStefano Zampini }
786c1f0200aSDmitry Karpeev 
7874d615ea0SBarry Smith /*@
788811af0c4SBarry Smith   PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
789811af0c4SBarry Smith   changes a second `PetscMPIInt` array to match the sorted first array.
7904d615ea0SBarry Smith 
7914d615ea0SBarry Smith   Not Collective
7924d615ea0SBarry Smith 
7934d615ea0SBarry Smith   Input Parameters:
7944d615ea0SBarry Smith + n - number of values
79559e16d97SJunchao Zhang . X - array of integers
79659e16d97SJunchao Zhang - Y - second array of integers
7974d615ea0SBarry Smith 
7984d615ea0SBarry Smith   Level: intermediate
7994d615ea0SBarry Smith 
800db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
8014d615ea0SBarry Smith @*/
802d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[])
803d71ae5a4SJacob Faibussowitsch {
80459e16d97SJunchao Zhang   PetscMPIInt pivot, t1, t2;
8054d615ea0SBarry Smith 
8064d615ea0SBarry Smith   PetscFunctionBegin;
8075f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
8083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8095569a785SJunchao Zhang }
8105569a785SJunchao Zhang 
8115569a785SJunchao Zhang /*@
812811af0c4SBarry Smith   PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
813811af0c4SBarry Smith   changes a second array of `PetscInt` to match the sorted first array.
8145569a785SJunchao Zhang 
8155569a785SJunchao Zhang   Not Collective
8165569a785SJunchao Zhang 
8175569a785SJunchao Zhang   Input Parameters:
8185569a785SJunchao Zhang + n - number of values
81959e16d97SJunchao Zhang . X - array of MPI integers
82059e16d97SJunchao Zhang - Y - second array of Petsc integers
8215569a785SJunchao Zhang 
8225569a785SJunchao Zhang   Level: intermediate
8235569a785SJunchao Zhang 
824811af0c4SBarry Smith   Note:
825811af0c4SBarry Smith   This routine is useful when one needs to sort MPI ranks with other integer arrays.
8265569a785SJunchao Zhang 
827db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
8285569a785SJunchao Zhang @*/
829d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[])
830d71ae5a4SJacob Faibussowitsch {
83159e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
8325569a785SJunchao Zhang   PetscInt    t2;
8335569a785SJunchao Zhang 
8345569a785SJunchao Zhang   PetscFunctionBegin;
8355f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
8363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
837e5c89e4eSSatish Balay }
838e5c89e4eSSatish Balay 
839e5c89e4eSSatish Balay /*@
840811af0c4SBarry Smith   PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
841811af0c4SBarry Smith   changes a second `PetscScalar` array to match the sorted first array.
842e5c89e4eSSatish Balay 
843e5c89e4eSSatish Balay   Not Collective
844e5c89e4eSSatish Balay 
845e5c89e4eSSatish Balay   Input Parameters:
846e5c89e4eSSatish Balay + n - number of values
84759e16d97SJunchao Zhang . X - array of integers
84859e16d97SJunchao Zhang - Y - second array of scalars
849e5c89e4eSSatish Balay 
850e5c89e4eSSatish Balay   Level: intermediate
851e5c89e4eSSatish Balay 
852db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
853e5c89e4eSSatish Balay @*/
854d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[])
855d71ae5a4SJacob Faibussowitsch {
85659e16d97SJunchao Zhang   PetscInt    pivot, t1;
85759e16d97SJunchao Zhang   PetscScalar t2;
858e5c89e4eSSatish Balay 
859e5c89e4eSSatish Balay   PetscFunctionBegin;
8605f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
8613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86217df18f2SToby Isaac }
86317df18f2SToby Isaac 
8646c2863d0SBarry Smith /*@C
865811af0c4SBarry Smith   PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
86617df18f2SToby Isaac   changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
86717df18f2SToby Isaac   provide workspace (the size of an element in the data array) to use when sorting.
86817df18f2SToby Isaac 
86917df18f2SToby Isaac   Not Collective
87017df18f2SToby Isaac 
87117df18f2SToby Isaac   Input Parameters:
87217df18f2SToby Isaac + n    - number of values
87359e16d97SJunchao Zhang . X    - array of integers
87459e16d97SJunchao Zhang . Y    - second array of data
87517df18f2SToby Isaac . size - sizeof elements in the data array in bytes
87659e16d97SJunchao Zhang - t2   - workspace of "size" bytes used when sorting
87717df18f2SToby Isaac 
87817df18f2SToby Isaac   Level: intermediate
87917df18f2SToby Isaac 
880db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
88117df18f2SToby Isaac @*/
882d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2)
883d71ae5a4SJacob Faibussowitsch {
88459e16d97SJunchao Zhang   char    *YY = (char *)Y;
8855f80ce2aSJacob Faibussowitsch   PetscInt t1, pivot, hi = n - 1;
88617df18f2SToby Isaac 
88717df18f2SToby Isaac   PetscFunctionBegin;
88817df18f2SToby Isaac   if (n < 8) {
8895f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < n; i++) {
89059e16d97SJunchao Zhang       pivot = X[i];
8915f80ce2aSJacob Faibussowitsch       for (PetscInt j = i + 1; j < n; j++) {
89259e16d97SJunchao Zhang         if (pivot > X[j]) {
89359e16d97SJunchao Zhang           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
89459e16d97SJunchao Zhang           pivot = X[i];
89517df18f2SToby Isaac         }
89617df18f2SToby Isaac       }
89717df18f2SToby Isaac     }
89817df18f2SToby Isaac   } else {
89959e16d97SJunchao Zhang     /* Two way partition */
9005f80ce2aSJacob Faibussowitsch     PetscInt l = 0, r = hi;
9015f80ce2aSJacob Faibussowitsch 
9025f80ce2aSJacob Faibussowitsch     pivot = X[MEDIAN(X, hi)];
90359e16d97SJunchao Zhang     while (1) {
90459e16d97SJunchao Zhang       while (X[l] < pivot) l++;
90559e16d97SJunchao Zhang       while (X[r] > pivot) r--;
9069371c9d4SSatish Balay       if (l >= r) {
9079371c9d4SSatish Balay         r++;
9089371c9d4SSatish Balay         break;
9099371c9d4SSatish Balay       }
91059e16d97SJunchao Zhang       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
91159e16d97SJunchao Zhang       l++;
91259e16d97SJunchao Zhang       r--;
91359e16d97SJunchao Zhang     }
9149566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
9159566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
91617df18f2SToby Isaac   }
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91817df18f2SToby Isaac }
91917df18f2SToby Isaac 
92021e72a00SBarry Smith /*@
921811af0c4SBarry Smith   PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
92221e72a00SBarry Smith 
92321e72a00SBarry Smith   Not Collective
92421e72a00SBarry Smith 
92521e72a00SBarry Smith   Input Parameters:
92621e72a00SBarry Smith + an - number of values in the first array
92721e72a00SBarry Smith . aI - first sorted array of integers
92821e72a00SBarry Smith . bn - number of values in the second array
92921e72a00SBarry Smith - bI - second array of integers
93021e72a00SBarry Smith 
93121e72a00SBarry Smith   Output Parameters:
93221e72a00SBarry Smith + n - number of values in the merged array
9336c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided
93421e72a00SBarry Smith 
93521e72a00SBarry Smith   Level: intermediate
93621e72a00SBarry Smith 
937db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
93821e72a00SBarry Smith @*/
939d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
940d71ae5a4SJacob Faibussowitsch {
94121e72a00SBarry Smith   PetscInt *L_ = *L, ak, bk, k;
94221e72a00SBarry Smith 
943362febeeSStefano Zampini   PetscFunctionBegin;
94421e72a00SBarry Smith   if (!L_) {
9459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(an + bn, L));
94621e72a00SBarry Smith     L_ = *L;
94721e72a00SBarry Smith   }
94821e72a00SBarry Smith   k = ak = bk = 0;
94921e72a00SBarry Smith   while (ak < an && bk < bn) {
95021e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
95121e72a00SBarry Smith       L_[k] = aI[ak];
95221e72a00SBarry Smith       ++ak;
95321e72a00SBarry Smith       ++bk;
95421e72a00SBarry Smith       ++k;
95521e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
95621e72a00SBarry Smith       L_[k] = aI[ak];
95721e72a00SBarry Smith       ++ak;
95821e72a00SBarry Smith       ++k;
95921e72a00SBarry Smith     } else {
96021e72a00SBarry Smith       L_[k] = bI[bk];
96121e72a00SBarry Smith       ++bk;
96221e72a00SBarry Smith       ++k;
96321e72a00SBarry Smith     }
96421e72a00SBarry Smith   }
96521e72a00SBarry Smith   if (ak < an) {
9669566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
96721e72a00SBarry Smith     k += (an - ak);
96821e72a00SBarry Smith   }
96921e72a00SBarry Smith   if (bk < bn) {
9709566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
97121e72a00SBarry Smith     k += (bn - bk);
97221e72a00SBarry Smith   }
97321e72a00SBarry Smith   *n = k;
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97521e72a00SBarry Smith }
976b4301105SBarry Smith 
977c1f0200aSDmitry Karpeev /*@
978811af0c4SBarry Smith   PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
979c1f0200aSDmitry Karpeev   The additional arrays are the same length as sorted arrays and are merged
980c1f0200aSDmitry Karpeev   in the order determined by the merging of the sorted pair.
981c1f0200aSDmitry Karpeev 
982c1f0200aSDmitry Karpeev   Not Collective
983c1f0200aSDmitry Karpeev 
984c1f0200aSDmitry Karpeev   Input Parameters:
985c1f0200aSDmitry Karpeev + an - number of values in the first array
986c1f0200aSDmitry Karpeev . aI - first sorted array of integers
987c1f0200aSDmitry Karpeev . aJ - first additional array of integers
988c1f0200aSDmitry Karpeev . bn - number of values in the second array
989c1f0200aSDmitry Karpeev . bI - second array of integers
990c1f0200aSDmitry Karpeev - bJ - second additional of integers
991c1f0200aSDmitry Karpeev 
992c1f0200aSDmitry Karpeev   Output Parameters:
993c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn)
99414c49006SJed Brown . L - merged sorted array
995c1f0200aSDmitry Karpeev - J - merged additional array
996c1f0200aSDmitry Karpeev 
997811af0c4SBarry Smith   Note:
998da81f932SPierre 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
999811af0c4SBarry Smith 
1000c1f0200aSDmitry Karpeev   Level: intermediate
1001c1f0200aSDmitry Karpeev 
1002db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1003c1f0200aSDmitry Karpeev @*/
1004d71ae5a4SJacob 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)
1005d71ae5a4SJacob Faibussowitsch {
100670d8d27cSBarry Smith   PetscInt n_, *L_, *J_, ak, bk, k;
1007c1f0200aSDmitry Karpeev 
100814c49006SJed Brown   PetscFunctionBegin;
10094f572ea9SToby Isaac   PetscAssertPointer(L, 8);
10104f572ea9SToby Isaac   PetscAssertPointer(J, 9);
1011c1f0200aSDmitry Karpeev   n_ = an + bn;
1012c1f0200aSDmitry Karpeev   *n = n_;
101348a46eb9SPierre Jolivet   if (!*L) PetscCall(PetscMalloc1(n_, L));
1014d7aa01a8SHong Zhang   L_ = *L;
101548a46eb9SPierre Jolivet   if (!*J) PetscCall(PetscMalloc1(n_, J));
1016c1f0200aSDmitry Karpeev   J_ = *J;
1017c1f0200aSDmitry Karpeev   k = ak = bk = 0;
1018c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
1019c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
1020d7aa01a8SHong Zhang       L_[k] = aI[ak];
1021c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
1022c1f0200aSDmitry Karpeev       ++ak;
1023c1f0200aSDmitry Karpeev       ++k;
1024a297a907SKarl Rupp     } else {
1025d7aa01a8SHong Zhang       L_[k] = bI[bk];
1026c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
1027c1f0200aSDmitry Karpeev       ++bk;
1028c1f0200aSDmitry Karpeev       ++k;
1029c1f0200aSDmitry Karpeev     }
1030c1f0200aSDmitry Karpeev   }
1031c1f0200aSDmitry Karpeev   if (ak < an) {
10329566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
10339566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1034c1f0200aSDmitry Karpeev     k += (an - ak);
1035c1f0200aSDmitry Karpeev   }
1036c1f0200aSDmitry Karpeev   if (bk < bn) {
10379566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
10389566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1039c1f0200aSDmitry Karpeev   }
10403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1041c1f0200aSDmitry Karpeev }
1042c1f0200aSDmitry Karpeev 
1043e498c390SJed Brown /*@
1044811af0c4SBarry Smith   PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
1045e498c390SJed Brown 
1046e498c390SJed Brown   Not Collective
1047e498c390SJed Brown 
1048e498c390SJed Brown   Input Parameters:
1049e498c390SJed Brown + an - number of values in the first array
1050e498c390SJed Brown . aI - first sorted array of integers
1051e498c390SJed Brown . bn - number of values in the second array
1052e498c390SJed Brown - bI - second array of integers
1053e498c390SJed Brown 
1054e498c390SJed Brown   Output Parameters:
1055e498c390SJed Brown + n - number of values in the merged array (<= an + bn)
1056e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed
1057e498c390SJed Brown 
1058e498c390SJed Brown   Level: intermediate
1059e498c390SJed Brown 
1060db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1061e498c390SJed Brown @*/
1062d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1063d71ae5a4SJacob Faibussowitsch {
1064e498c390SJed Brown   PetscInt ai, bi, k;
1065e498c390SJed Brown 
1066e498c390SJed Brown   PetscFunctionBegin;
10679566063dSJacob Faibussowitsch   if (!*L) PetscCall(PetscMalloc1((an + bn), L));
1068e498c390SJed Brown   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1069e498c390SJed Brown     PetscInt t = -1;
1070e498c390SJed Brown     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
1071fbccb6d4SPierre Jolivet     for (; bi < bn && bI[bi] == t; bi++);
1072e498c390SJed Brown     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
1073fbccb6d4SPierre Jolivet     for (; ai < an && aI[ai] == t; ai++);
1074e498c390SJed Brown   }
1075e498c390SJed Brown   *n = k;
10763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1077e498c390SJed Brown }
1078e5c89e4eSSatish Balay 
10796c2863d0SBarry Smith /*@C
108042eaa7f3SBarry Smith   PetscProcessTree - Prepares tree data to be displayed graphically
108142eaa7f3SBarry Smith 
108242eaa7f3SBarry Smith   Not Collective
108342eaa7f3SBarry Smith 
108442eaa7f3SBarry Smith   Input Parameters:
108542eaa7f3SBarry Smith + n        - number of values
108642eaa7f3SBarry Smith . mask     - indicates those entries in the tree, location 0 is always masked
108742eaa7f3SBarry Smith - parentid - indicates the parent of each entry
108842eaa7f3SBarry Smith 
108942eaa7f3SBarry Smith   Output Parameters:
109042eaa7f3SBarry Smith + Nlevels   - the number of levels
109142eaa7f3SBarry Smith . Level     - for each node tells its level
1092aec76313SJacob Faibussowitsch . Levelcnt  - the number of nodes on each level
109342eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc
109442eaa7f3SBarry Smith - Column    - for each id tells its column index
109542eaa7f3SBarry Smith 
10966c2863d0SBarry Smith   Level: developer
109742eaa7f3SBarry Smith 
1098811af0c4SBarry Smith   Note:
109995452b02SPatrick Sanan   This code is not currently used
110042eaa7f3SBarry Smith 
1101db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
110242eaa7f3SBarry Smith @*/
1103d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1104d71ae5a4SJacob Faibussowitsch {
110542eaa7f3SBarry Smith   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1106ace3abfcSBarry Smith   PetscBool done = PETSC_FALSE;
110742eaa7f3SBarry Smith 
110842eaa7f3SBarry Smith   PetscFunctionBegin;
110908401ef6SPierre Jolivet   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
111042eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
111142eaa7f3SBarry Smith     if (mask[i]) continue;
111208401ef6SPierre Jolivet     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1113cc73adaaSBarry Smith     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
111442eaa7f3SBarry Smith   }
111542eaa7f3SBarry Smith 
111642eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
111742eaa7f3SBarry Smith     if (!mask[i]) nmask++;
111842eaa7f3SBarry Smith   }
111942eaa7f3SBarry Smith 
112042eaa7f3SBarry Smith   /* determine the level in the tree of each node */
11219566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &level));
1122a297a907SKarl Rupp 
112342eaa7f3SBarry Smith   level[0] = 1;
112442eaa7f3SBarry Smith   while (!done) {
112542eaa7f3SBarry Smith     done = PETSC_TRUE;
112642eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
112742eaa7f3SBarry Smith       if (mask[i]) continue;
112842eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
112942eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
113042eaa7f3SBarry Smith     }
113142eaa7f3SBarry Smith   }
113242eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
113342eaa7f3SBarry Smith     level[i]--;
113442eaa7f3SBarry Smith     nlevels = PetscMax(nlevels, level[i]);
113542eaa7f3SBarry Smith   }
113642eaa7f3SBarry Smith 
113742eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
11389566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nlevels, &levelcnt));
113942eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
114042eaa7f3SBarry Smith     if (mask[i]) continue;
114142eaa7f3SBarry Smith     levelcnt[level[i] - 1]++;
114242eaa7f3SBarry Smith   }
1143a297a907SKarl Rupp   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
114442eaa7f3SBarry Smith 
114542eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
11469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
11479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmask, &idbylevel));
114842eaa7f3SBarry Smith   for (j = 1; j <= nlevels; j++) {
114942eaa7f3SBarry Smith     cnt = 0;
115042eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
115142eaa7f3SBarry Smith       if (mask[i]) continue;
115242eaa7f3SBarry Smith       if (level[i] != j) continue;
115342eaa7f3SBarry Smith       workid[cnt]         = i;
115442eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
115542eaa7f3SBarry Smith     }
115642eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
115742eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
11589566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
115942eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
116042eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
11619566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
116242eaa7f3SBarry Smith     tcnt += cnt;
116342eaa7f3SBarry Smith   }
116408401ef6SPierre Jolivet   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
11659566063dSJacob Faibussowitsch   PetscCall(PetscFree2(workid, workparentid));
116642eaa7f3SBarry Smith 
116742eaa7f3SBarry Smith   /* for each node list its column */
11689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &column));
116942eaa7f3SBarry Smith   cnt = 0;
117042eaa7f3SBarry Smith   for (j = 0; j < nlevels; j++) {
1171ad540459SPierre Jolivet     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
117242eaa7f3SBarry Smith   }
117342eaa7f3SBarry Smith 
117442eaa7f3SBarry Smith   *Nlevels   = nlevels;
117542eaa7f3SBarry Smith   *Level     = level;
117642eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
117742eaa7f3SBarry Smith   *Idbylevel = idbylevel;
117842eaa7f3SBarry Smith   *Column    = column;
11793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118042eaa7f3SBarry Smith }
1181ce605777SToby Isaac 
1182ce605777SToby Isaac /*@
1183811af0c4SBarry Smith   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1184ce605777SToby Isaac 
1185ce605777SToby Isaac   Collective
1186ce605777SToby Isaac 
1187ce605777SToby Isaac   Input Parameters:
1188ce605777SToby Isaac + comm - the MPI communicator
1189ce605777SToby Isaac . n    - the local number of integers
1190ce605777SToby Isaac - keys - the local array of integers
1191ce605777SToby Isaac 
1192ce605777SToby Isaac   Output Parameters:
1193ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1194ce605777SToby Isaac 
1195ce605777SToby Isaac   Level: developer
1196ce605777SToby Isaac 
1197db781477SPatrick Sanan .seealso: `PetscParallelSortInt()`
1198ce605777SToby Isaac @*/
1199d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1200d71ae5a4SJacob Faibussowitsch {
1201ce605777SToby Isaac   PetscBool   sorted;
1202ce605777SToby Isaac   PetscInt    i, min, max, prevmax;
12032a1da528SToby Isaac   PetscMPIInt rank;
1204ce605777SToby Isaac 
1205ce605777SToby Isaac   PetscFunctionBegin;
1206ce605777SToby Isaac   sorted = PETSC_TRUE;
1207ce605777SToby Isaac   min    = PETSC_MAX_INT;
1208ce605777SToby Isaac   max    = PETSC_MIN_INT;
1209ce605777SToby Isaac   if (n) {
1210ce605777SToby Isaac     min = keys[0];
1211ce605777SToby Isaac     max = keys[0];
1212ce605777SToby Isaac   }
1213ce605777SToby Isaac   for (i = 1; i < n; i++) {
1214ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1215ce605777SToby Isaac     min = PetscMin(min, keys[i]);
1216ce605777SToby Isaac     max = PetscMax(max, keys[i]);
1217ce605777SToby Isaac   }
1218ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1219ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
12209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
12219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1222dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1223ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1224712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
12253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1226ce605777SToby Isaac }
1227