xref: /petsc/src/sys/utils/sorti.c (revision 54c0599748d6bec40ddc3bedea0515243cbbf0fb)
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
2493334a463SJames Wright - X - array of `PetscInt`
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 @*/
2586497c311SBarry Smith PetscErrorCode PetscSortedInt(PetscCount 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
2743334a463SJames Wright - X - array of `PetscInt64`
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 @*/
2836497c311SBarry Smith PetscErrorCode PetscSortedInt64(PetscCount 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
2993334a463SJames Wright - X - array of `PetscInt`
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 @*/
3106497c311SBarry Smith PetscErrorCode PetscSortInt(PetscCount n, PetscInt X[])
311d71ae5a4SJacob Faibussowitsch {
3126497c311SBarry Smith   PetscInt pivot, t1, N;
313e5c89e4eSSatish Balay 
314e5c89e4eSSatish Balay   PetscFunctionBegin;
3154f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
3166497c311SBarry Smith   PetscCall(PetscIntCast(n, &N));
3176497c311SBarry Smith   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
3283334a463SJames Wright - X - array of `PetscInt64`
32962e5d2d2SJDBetteridge 
33062e5d2d2SJDBetteridge   Notes:
3317a43aa34SPierre Jolivet   This function sorts `PetscInt64`s assumed to be in completely random order
33262e5d2d2SJDBetteridge 
33362e5d2d2SJDBetteridge   Level: intermediate
33462e5d2d2SJDBetteridge 
33562e5d2d2SJDBetteridge .seealso: `PetscSortInt()`
33662e5d2d2SJDBetteridge @*/
3376497c311SBarry Smith PetscErrorCode PetscSortInt64(PetscCount n, PetscInt64 X[])
33862e5d2d2SJDBetteridge {
33962e5d2d2SJDBetteridge   PetscCount pivot, t1;
34062e5d2d2SJDBetteridge 
34162e5d2d2SJDBetteridge   PetscFunctionBegin;
3424f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
34362e5d2d2SJDBetteridge   QuickSort1(PetscSortInt64, X, n, pivot, t1);
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34562e5d2d2SJDBetteridge }
34662e5d2d2SJDBetteridge 
34762e5d2d2SJDBetteridge /*@
3483334a463SJames Wright   PetscSortCount - Sorts an array of `PetscCount` in place in increasing order.
34962e5d2d2SJDBetteridge 
35062e5d2d2SJDBetteridge   Not Collective
35162e5d2d2SJDBetteridge 
35262e5d2d2SJDBetteridge   Input Parameters:
35362e5d2d2SJDBetteridge + n - number of values
3543334a463SJames Wright - X - array of `PetscCount`
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 @*/
3636497c311SBarry Smith PetscErrorCode PetscSortCount(PetscCount n, PetscCount X[])
36462e5d2d2SJDBetteridge {
36562e5d2d2SJDBetteridge   PetscCount pivot, t1;
36662e5d2d2SJDBetteridge 
36762e5d2d2SJDBetteridge   PetscFunctionBegin;
3684f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
36962e5d2d2SJDBetteridge   QuickSort1(PetscSortCount, X, n, pivot, t1);
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37162e5d2d2SJDBetteridge }
37262e5d2d2SJDBetteridge 
37362e5d2d2SJDBetteridge /*@
3743334a463SJames Wright   PetscSortReverseInt - Sorts an array of `PetscInt` in place in decreasing order.
375ce605777SToby Isaac 
376ce605777SToby Isaac   Not Collective
377ce605777SToby Isaac 
378ce605777SToby Isaac   Input Parameters:
379ce605777SToby Isaac + n - number of values
3803334a463SJames Wright - X - array of `PetscInt`
381ce605777SToby Isaac 
382ce605777SToby Isaac   Level: intermediate
383ce605777SToby Isaac 
384db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
385ce605777SToby Isaac @*/
3866497c311SBarry Smith PetscErrorCode PetscSortReverseInt(PetscCount n, PetscInt X[])
387d71ae5a4SJacob Faibussowitsch {
388ce605777SToby Isaac   PetscInt pivot, t1;
389ce605777SToby Isaac 
390ce605777SToby Isaac   PetscFunctionBegin;
3914f572ea9SToby Isaac   if (n) PetscAssertPointer(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
4033334a463SJames Wright - X - sorted array of `PetscInt`
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;
4174f572ea9SToby Isaac   PetscAssertPointer(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
4363334a463SJames Wright - X - sorted array of `PetscInt`
437b6a18d21SVaclav Hapla 
438b6a18d21SVaclav Hapla   Output Parameter:
4393334a463SJames Wright . flg - True if the array has duplications, otherwise false
440b6a18d21SVaclav Hapla 
441b6a18d21SVaclav Hapla   Level: intermediate
442b6a18d21SVaclav Hapla 
443db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
444b6a18d21SVaclav Hapla @*/
4456497c311SBarry Smith PetscErrorCode PetscSortedCheckDupsInt(PetscCount 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
4683334a463SJames Wright - X - array of `PetscInt`
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;
4804f572ea9SToby Isaac   PetscAssertPointer(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 /*@
4873334a463SJames Wright   PetscFindInt - Finds the location of a `PetscInt` key in a sorted array of `PetscInt`
48860e03357SMatthew G Knepley 
48960e03357SMatthew G Knepley   Not Collective
49060e03357SMatthew G Knepley 
49160e03357SMatthew G Knepley   Input Parameters:
4923334a463SJames Wright + key - the `PetscInt` key to locate
493d9f1162dSJed Brown . n   - number of values in the array
4943334a463SJames Wright - X   - array of `PetscInt`
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 @*/
5036497c311SBarry Smith PetscErrorCode PetscFindInt(PetscInt key, PetscCount n, const PetscInt X[], PetscInt *loc)
504d71ae5a4SJacob Faibussowitsch {
5056497c311SBarry Smith   PetscInt lo = 0, hi;
50660e03357SMatthew G Knepley 
50760e03357SMatthew G Knepley   PetscFunctionBegin;
5084f572ea9SToby Isaac   PetscAssertPointer(loc, 4);
5099371c9d4SSatish Balay   if (!n) {
5109371c9d4SSatish Balay     *loc = -1;
5113ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5129371c9d4SSatish Balay   }
5134f572ea9SToby Isaac   PetscAssertPointer(X, 3);
5146497c311SBarry Smith   PetscCall(PetscIntCast(n, &hi));
515d9f1162dSJed Brown   while (hi - lo > 1) {
516d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo) / 2;
51742f965a0SMatthew 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]);
51859e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
519d9f1162dSJed Brown     else lo = mid;
52060e03357SMatthew G Knepley   }
52159e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52360e03357SMatthew G Knepley }
52460e03357SMatthew G Knepley 
525d2aeb606SJed Brown /*@
526811af0c4SBarry Smith   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
527f1cab4e1SJunchao Zhang 
528f1cab4e1SJunchao Zhang   Not Collective
529f1cab4e1SJunchao Zhang 
530f1cab4e1SJunchao Zhang   Input Parameters:
531f1cab4e1SJunchao Zhang + n - number of values in the array
5323334a463SJames Wright - X - array of `PetscInt`
533f1cab4e1SJunchao Zhang 
534f1cab4e1SJunchao Zhang   Output Parameter:
535f1cab4e1SJunchao Zhang . dups - True if the array has dups, otherwise false
536f1cab4e1SJunchao Zhang 
537f1cab4e1SJunchao Zhang   Level: intermediate
538f1cab4e1SJunchao Zhang 
539db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
540f1cab4e1SJunchao Zhang @*/
541d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
542d71ae5a4SJacob Faibussowitsch {
543f1cab4e1SJunchao Zhang   PetscInt   i;
544f1cab4e1SJunchao Zhang   PetscHSetI ht;
545f1cab4e1SJunchao Zhang   PetscBool  missing;
546f1cab4e1SJunchao Zhang 
547f1cab4e1SJunchao Zhang   PetscFunctionBegin;
5484f572ea9SToby Isaac   if (n) PetscAssertPointer(X, 2);
5494f572ea9SToby Isaac   PetscAssertPointer(dups, 3);
550f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
551f1cab4e1SJunchao Zhang   if (n > 1) {
5529566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&ht));
5539566063dSJacob Faibussowitsch     PetscCall(PetscHSetIResize(ht, n));
554f1cab4e1SJunchao Zhang     for (i = 0; i < n; i++) {
5559566063dSJacob Faibussowitsch       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
5569371c9d4SSatish Balay       if (!missing) {
5579371c9d4SSatish Balay         *dups = PETSC_TRUE;
5589371c9d4SSatish Balay         break;
5599371c9d4SSatish Balay       }
560f1cab4e1SJunchao Zhang     }
5619566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&ht));
562f1cab4e1SJunchao Zhang   }
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
564f1cab4e1SJunchao Zhang }
565f1cab4e1SJunchao Zhang 
566f1cab4e1SJunchao Zhang /*@
567811af0c4SBarry Smith   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
568d2aeb606SJed Brown 
569d2aeb606SJed Brown   Not Collective
570d2aeb606SJed Brown 
571d2aeb606SJed Brown   Input Parameters:
572d2aeb606SJed Brown + key - the integer to locate
573d2aeb606SJed Brown . n   - number of values in the array
5743334a463SJames Wright - X   - array of `PetscMPIInt`
575d2aeb606SJed Brown 
576d2aeb606SJed Brown   Output Parameter:
577d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
578d2aeb606SJed Brown 
579d2aeb606SJed Brown   Level: intermediate
580d2aeb606SJed Brown 
581db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
582d2aeb606SJed Brown @*/
5836497c311SBarry Smith PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscCount n, const PetscMPIInt X[], PetscInt *loc)
584d71ae5a4SJacob Faibussowitsch {
5856497c311SBarry Smith   PetscCount lo = 0, hi = n;
586d2aeb606SJed Brown 
587d2aeb606SJed Brown   PetscFunctionBegin;
5884f572ea9SToby Isaac   PetscAssertPointer(loc, 4);
5899371c9d4SSatish Balay   if (!n) {
5909371c9d4SSatish Balay     *loc = -1;
5913ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5929371c9d4SSatish Balay   }
5934f572ea9SToby Isaac   PetscAssertPointer(X, 3);
594d2aeb606SJed Brown   while (hi - lo > 1) {
5956497c311SBarry Smith     PetscCount mid = lo + (hi - lo) / 2;
59642f965a0SMatthew 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]);
59759e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
598d2aeb606SJed Brown     else lo = mid;
599d2aeb606SJed Brown   }
6006497c311SBarry Smith   PetscCall(PetscIntCast(key == X[lo] ? lo : -(lo + (key > X[lo]) + 1), loc));
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
602e5c89e4eSSatish Balay }
603e5c89e4eSSatish Balay 
604e5c89e4eSSatish Balay /*@
605811af0c4SBarry Smith   PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
606811af0c4SBarry Smith   changes a second array of `PetscInt` to match the sorted first array.
607e5c89e4eSSatish Balay 
608e5c89e4eSSatish Balay   Not Collective
609e5c89e4eSSatish Balay 
610e5c89e4eSSatish Balay   Input Parameters:
611e5c89e4eSSatish Balay + n - number of values
6123334a463SJames Wright . X - array of `PetscInt`
6133334a463SJames Wright - Y - second array of `PetscInt`
614e5c89e4eSSatish Balay 
615e5c89e4eSSatish Balay   Level: intermediate
616e5c89e4eSSatish Balay 
617db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
618e5c89e4eSSatish Balay @*/
6196497c311SBarry Smith PetscErrorCode PetscSortIntWithArray(PetscCount n, PetscInt X[], PetscInt Y[])
620d71ae5a4SJacob Faibussowitsch {
62159e16d97SJunchao Zhang   PetscInt pivot, t1, t2;
622e5c89e4eSSatish Balay 
623e5c89e4eSSatish Balay   PetscFunctionBegin;
6245f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
626c1f0200aSDmitry Karpeev }
627c1f0200aSDmitry Karpeev 
628c1f0200aSDmitry Karpeev /*@
629811af0c4SBarry Smith   PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
630811af0c4SBarry Smith   changes a pair of `PetscInt` arrays to match the sorted first array.
631c1f0200aSDmitry Karpeev 
632c1f0200aSDmitry Karpeev   Not Collective
633c1f0200aSDmitry Karpeev 
634c1f0200aSDmitry Karpeev   Input Parameters:
635c1f0200aSDmitry Karpeev + n - number of values
6363334a463SJames Wright . X - array of `PestcInt`
6373334a463SJames Wright . Y - second array of `PestcInt` (first array of the pair)
6383334a463SJames Wright - Z - third array of `PestcInt` (second array of the pair)
639c1f0200aSDmitry Karpeev 
640c1f0200aSDmitry Karpeev   Level: intermediate
641c1f0200aSDmitry Karpeev 
642db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
643c1f0200aSDmitry Karpeev @*/
6446497c311SBarry Smith PetscErrorCode PetscSortIntWithArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscInt Z[])
645d71ae5a4SJacob Faibussowitsch {
64659e16d97SJunchao Zhang   PetscInt pivot, t1, t2, t3;
647c1f0200aSDmitry Karpeev 
648c1f0200aSDmitry Karpeev   PetscFunctionBegin;
6495f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
651c1f0200aSDmitry Karpeev }
652c1f0200aSDmitry Karpeev 
65317d7d925SStefano Zampini /*@
6546497c311SBarry Smith   PetscSortIntWithMPIIntArray - Sorts an array of `PetscInt` in place in increasing order;
6556497c311SBarry Smith   changes a second array of `PetscMPI` to match the sorted first array.
6566497c311SBarry Smith 
6576497c311SBarry Smith   Not Collective
6586497c311SBarry Smith 
6596497c311SBarry Smith   Input Parameters:
6606497c311SBarry Smith + n - number of values
6613334a463SJames Wright . X - array of `PetscInt`
6623334a463SJames Wright - Y - second array of `PetscMPIInt`
6636497c311SBarry Smith 
6646497c311SBarry Smith   Level: intermediate
6656497c311SBarry Smith 
666*54c05997SPierre Jolivet .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
6676497c311SBarry Smith @*/
6686497c311SBarry Smith PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount n, PetscInt X[], PetscMPIInt Y[])
6696497c311SBarry Smith {
6706497c311SBarry Smith   PetscInt    pivot, t1;
6716497c311SBarry Smith   PetscMPIInt t2;
6726497c311SBarry Smith 
6736497c311SBarry Smith   PetscFunctionBegin;
6746497c311SBarry Smith   QuickSort2(PetscSortIntWithMPIIntArray, X, Y, n, pivot, t1, t2);
6756497c311SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6766497c311SBarry Smith }
6776497c311SBarry Smith 
6786497c311SBarry Smith /*@
679811af0c4SBarry Smith   PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
680811af0c4SBarry Smith   changes a second array of `PetscCount` to match the sorted first array.
681981bb840SJunchao Zhang 
682981bb840SJunchao Zhang   Not Collective
683981bb840SJunchao Zhang 
684981bb840SJunchao Zhang   Input Parameters:
685981bb840SJunchao Zhang + n - number of values
6863334a463SJames Wright . X - array of `PetscInt`
6873334a463SJames Wright - Y - second array of `PetscCount`
688981bb840SJunchao Zhang 
689981bb840SJunchao Zhang   Level: intermediate
690981bb840SJunchao Zhang 
691*54c05997SPierre Jolivet .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
692981bb840SJunchao Zhang @*/
693d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
694d71ae5a4SJacob Faibussowitsch {
695981bb840SJunchao Zhang   PetscInt   pivot, t1;
696981bb840SJunchao Zhang   PetscCount t2;
697981bb840SJunchao Zhang 
698981bb840SJunchao Zhang   PetscFunctionBegin;
6995f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
7003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
701981bb840SJunchao Zhang }
702981bb840SJunchao Zhang 
703981bb840SJunchao Zhang /*@
704811af0c4SBarry Smith   PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
705811af0c4SBarry Smith   changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
706981bb840SJunchao Zhang 
707981bb840SJunchao Zhang   Not Collective
708981bb840SJunchao Zhang 
709981bb840SJunchao Zhang   Input Parameters:
710981bb840SJunchao Zhang + n - number of values
7113334a463SJames Wright . X - array of `PetscInt`
7123334a463SJames Wright . Y - second array of `PetscInt` (first array of the pair)
7133334a463SJames Wright - Z - third array of `PetscCount` (second array of the pair)
714981bb840SJunchao Zhang 
715981bb840SJunchao Zhang   Level: intermediate
716981bb840SJunchao Zhang 
717811af0c4SBarry Smith   Note:
718981bb840SJunchao 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.
719981bb840SJunchao Zhang 
720*54c05997SPierre Jolivet .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
721981bb840SJunchao Zhang @*/
722d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
723d71ae5a4SJacob Faibussowitsch {
724981bb840SJunchao Zhang   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
725981bb840SJunchao Zhang   PetscCount t3;            /* temp for Z[] */
726981bb840SJunchao Zhang 
727981bb840SJunchao Zhang   PetscFunctionBegin;
7285f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
7293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
730981bb840SJunchao Zhang }
731981bb840SJunchao Zhang 
732981bb840SJunchao Zhang /*@
733811af0c4SBarry Smith   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
7346a7c706bSVaclav Hapla 
7356a7c706bSVaclav Hapla   Not Collective
7366a7c706bSVaclav Hapla 
7376a7c706bSVaclav Hapla   Input Parameters:
7386a7c706bSVaclav Hapla + n - number of values
7393334a463SJames Wright - X - array of `PetscMPIInt`
7406a7c706bSVaclav Hapla 
7412fe279fdSBarry Smith   Output Parameter:
7426a7c706bSVaclav Hapla . sorted - flag whether the array is sorted
7436a7c706bSVaclav Hapla 
7446a7c706bSVaclav Hapla   Level: intermediate
7456a7c706bSVaclav Hapla 
746db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
7476a7c706bSVaclav Hapla @*/
7486497c311SBarry Smith PetscErrorCode PetscSortedMPIInt(PetscCount n, const PetscMPIInt X[], PetscBool *sorted)
749d71ae5a4SJacob Faibussowitsch {
7506a7c706bSVaclav Hapla   PetscFunctionBegin;
7516a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
7523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7536a7c706bSVaclav Hapla }
7546a7c706bSVaclav Hapla 
7556a7c706bSVaclav Hapla /*@
756811af0c4SBarry Smith   PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
75717d7d925SStefano Zampini 
75817d7d925SStefano Zampini   Not Collective
75917d7d925SStefano Zampini 
76017d7d925SStefano Zampini   Input Parameters:
76117d7d925SStefano Zampini + n - number of values
7623334a463SJames Wright - X - array of `PetscMPIInt`
76317d7d925SStefano Zampini 
76417d7d925SStefano Zampini   Level: intermediate
76517d7d925SStefano Zampini 
766811af0c4SBarry Smith   Note:
767676f2a66SJacob Faibussowitsch   This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
768a5b23f4aSJose E. Roman   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
769676f2a66SJacob Faibussowitsch   code to see which routine is fastest.
770676f2a66SJacob Faibussowitsch 
771db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
77217d7d925SStefano Zampini @*/
7736497c311SBarry Smith PetscErrorCode PetscSortMPIInt(PetscCount n, PetscMPIInt X[])
774d71ae5a4SJacob Faibussowitsch {
77559e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
77617d7d925SStefano Zampini 
77717d7d925SStefano Zampini   PetscFunctionBegin;
7785f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
7793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78017d7d925SStefano Zampini }
78117d7d925SStefano Zampini 
78217d7d925SStefano Zampini /*@
783811af0c4SBarry Smith   PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
78417d7d925SStefano Zampini 
78517d7d925SStefano Zampini   Not Collective
78617d7d925SStefano Zampini 
78717d7d925SStefano Zampini   Input Parameters:
78817d7d925SStefano Zampini + n - number of values
7893334a463SJames Wright - X - array of `PetscMPIInt`
79017d7d925SStefano Zampini 
79117d7d925SStefano Zampini   Output Parameter:
79217d7d925SStefano Zampini . n - number of non-redundant values
79317d7d925SStefano Zampini 
79417d7d925SStefano Zampini   Level: intermediate
79517d7d925SStefano Zampini 
796db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
79717d7d925SStefano Zampini @*/
798d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
799d71ae5a4SJacob Faibussowitsch {
8005f80ce2aSJacob Faibussowitsch   PetscInt s = 0, N = *n, b = 0;
80117d7d925SStefano Zampini 
80217d7d925SStefano Zampini   PetscFunctionBegin;
8039566063dSJacob Faibussowitsch   PetscCall(PetscSortMPIInt(N, X));
8045f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < N - 1; i++) {
80559e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
8069371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
8079371c9d4SSatish Balay       b++;
808a297a907SKarl Rupp     } else s++;
80917d7d925SStefano Zampini   }
81017d7d925SStefano Zampini   *n = N - s;
8113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81217d7d925SStefano Zampini }
813c1f0200aSDmitry Karpeev 
8144d615ea0SBarry Smith /*@
815811af0c4SBarry Smith   PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
816811af0c4SBarry Smith   changes a second `PetscMPIInt` array to match the sorted first array.
8174d615ea0SBarry Smith 
8184d615ea0SBarry Smith   Not Collective
8194d615ea0SBarry Smith 
8204d615ea0SBarry Smith   Input Parameters:
8214d615ea0SBarry Smith + n - number of values
8223334a463SJames Wright . X - array of `PetscMPIInt`
8233334a463SJames Wright - Y - second array of `PetscMPIInt`
8244d615ea0SBarry Smith 
8254d615ea0SBarry Smith   Level: intermediate
8264d615ea0SBarry Smith 
827db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
8284d615ea0SBarry Smith @*/
8296497c311SBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscCount n, PetscMPIInt X[], PetscMPIInt Y[])
830d71ae5a4SJacob Faibussowitsch {
83159e16d97SJunchao Zhang   PetscMPIInt pivot, t1, t2;
8324d615ea0SBarry Smith 
8334d615ea0SBarry Smith   PetscFunctionBegin;
8345f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8365569a785SJunchao Zhang }
8375569a785SJunchao Zhang 
8385569a785SJunchao Zhang /*@
839811af0c4SBarry Smith   PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
840811af0c4SBarry Smith   changes a second array of `PetscInt` to match the sorted first array.
8415569a785SJunchao Zhang 
8425569a785SJunchao Zhang   Not Collective
8435569a785SJunchao Zhang 
8445569a785SJunchao Zhang   Input Parameters:
8455569a785SJunchao Zhang + n - number of values
8463334a463SJames Wright . X - array of `PetscMPIInt`
8473334a463SJames Wright - Y - second array of `PetscInt`
8485569a785SJunchao Zhang 
8495569a785SJunchao Zhang   Level: intermediate
8505569a785SJunchao Zhang 
851811af0c4SBarry Smith   Note:
852811af0c4SBarry Smith   This routine is useful when one needs to sort MPI ranks with other integer arrays.
8535569a785SJunchao Zhang 
854db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
8555569a785SJunchao Zhang @*/
8566497c311SBarry Smith PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount n, PetscMPIInt X[], PetscInt Y[])
857d71ae5a4SJacob Faibussowitsch {
85859e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
8595569a785SJunchao Zhang   PetscInt    t2;
8605569a785SJunchao Zhang 
8615569a785SJunchao Zhang   PetscFunctionBegin;
8625f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
864e5c89e4eSSatish Balay }
865e5c89e4eSSatish Balay 
866e5c89e4eSSatish Balay /*@
867811af0c4SBarry Smith   PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
868811af0c4SBarry Smith   changes a second `PetscScalar` array to match the sorted first array.
869e5c89e4eSSatish Balay 
870e5c89e4eSSatish Balay   Not Collective
871e5c89e4eSSatish Balay 
872e5c89e4eSSatish Balay   Input Parameters:
873e5c89e4eSSatish Balay + n - number of values
8743334a463SJames Wright . X - array of `PetscInt`
8753334a463SJames Wright - Y - second array of `PetscScalar`
876e5c89e4eSSatish Balay 
877e5c89e4eSSatish Balay   Level: intermediate
878e5c89e4eSSatish Balay 
879db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
880e5c89e4eSSatish Balay @*/
8816497c311SBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscCount n, PetscInt X[], PetscScalar Y[])
882d71ae5a4SJacob Faibussowitsch {
88359e16d97SJunchao Zhang   PetscInt    pivot, t1;
88459e16d97SJunchao Zhang   PetscScalar t2;
885e5c89e4eSSatish Balay 
886e5c89e4eSSatish Balay   PetscFunctionBegin;
8875f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
8883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88917df18f2SToby Isaac }
89017df18f2SToby Isaac 
8916c2863d0SBarry Smith /*@C
892811af0c4SBarry Smith   PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
89317df18f2SToby Isaac   changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
89417df18f2SToby Isaac   provide workspace (the size of an element in the data array) to use when sorting.
89517df18f2SToby Isaac 
896cc4c1da9SBarry Smith   Not Collective, No Fortran Support
89717df18f2SToby Isaac 
89817df18f2SToby Isaac   Input Parameters:
89917df18f2SToby Isaac + n    - number of values
9003334a463SJames Wright . X    - array of `PetscInt`
90159e16d97SJunchao Zhang . Y    - second array of data
90217df18f2SToby Isaac . size - sizeof elements in the data array in bytes
90359e16d97SJunchao Zhang - t2   - workspace of "size" bytes used when sorting
90417df18f2SToby Isaac 
90517df18f2SToby Isaac   Level: intermediate
90617df18f2SToby Isaac 
907db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
90817df18f2SToby Isaac @*/
9096497c311SBarry Smith PetscErrorCode PetscSortIntWithDataArray(PetscCount n, PetscInt X[], void *Y, size_t size, void *t2)
910d71ae5a4SJacob Faibussowitsch {
91159e16d97SJunchao Zhang   char      *YY = (char *)Y;
9126497c311SBarry Smith   PetscCount hi = n - 1;
9136497c311SBarry Smith   PetscInt   pivot, t1;
91417df18f2SToby Isaac 
91517df18f2SToby Isaac   PetscFunctionBegin;
91617df18f2SToby Isaac   if (n < 8) {
9176497c311SBarry Smith     for (PetscCount i = 0; i < n; i++) {
91859e16d97SJunchao Zhang       pivot = X[i];
9196497c311SBarry Smith       for (PetscCount j = i + 1; j < n; j++) {
92059e16d97SJunchao Zhang         if (pivot > X[j]) {
92159e16d97SJunchao Zhang           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
92259e16d97SJunchao Zhang           pivot = X[i];
92317df18f2SToby Isaac         }
92417df18f2SToby Isaac       }
92517df18f2SToby Isaac     }
92617df18f2SToby Isaac   } else {
92759e16d97SJunchao Zhang     /* Two way partition */
9286497c311SBarry Smith     PetscCount l = 0, r = hi;
9295f80ce2aSJacob Faibussowitsch 
9305f80ce2aSJacob Faibussowitsch     pivot = X[MEDIAN(X, hi)];
93159e16d97SJunchao Zhang     while (1) {
93259e16d97SJunchao Zhang       while (X[l] < pivot) l++;
93359e16d97SJunchao Zhang       while (X[r] > pivot) r--;
9349371c9d4SSatish Balay       if (l >= r) {
9359371c9d4SSatish Balay         r++;
9369371c9d4SSatish Balay         break;
9379371c9d4SSatish Balay       }
93859e16d97SJunchao Zhang       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
93959e16d97SJunchao Zhang       l++;
94059e16d97SJunchao Zhang       r--;
94159e16d97SJunchao Zhang     }
9429566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
9439566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
94417df18f2SToby Isaac   }
9453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
94617df18f2SToby Isaac }
94717df18f2SToby Isaac 
94821e72a00SBarry Smith /*@
949811af0c4SBarry Smith   PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
95021e72a00SBarry Smith 
95121e72a00SBarry Smith   Not Collective
95221e72a00SBarry Smith 
95321e72a00SBarry Smith   Input Parameters:
95421e72a00SBarry Smith + an - number of values in the first array
9553334a463SJames Wright . aI - first sorted array of `PetscInt`
95621e72a00SBarry Smith . bn - number of values in the second array
9573334a463SJames Wright - bI - second array of `PetscInt`
95821e72a00SBarry Smith 
95921e72a00SBarry Smith   Output Parameters:
96021e72a00SBarry Smith + n - number of values in the merged array
9616c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided
96221e72a00SBarry Smith 
96321e72a00SBarry Smith   Level: intermediate
96421e72a00SBarry Smith 
965db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
96621e72a00SBarry Smith @*/
967d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
968d71ae5a4SJacob Faibussowitsch {
96921e72a00SBarry Smith   PetscInt *L_ = *L, ak, bk, k;
97021e72a00SBarry Smith 
971362febeeSStefano Zampini   PetscFunctionBegin;
97221e72a00SBarry Smith   if (!L_) {
9739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(an + bn, L));
97421e72a00SBarry Smith     L_ = *L;
97521e72a00SBarry Smith   }
97621e72a00SBarry Smith   k = ak = bk = 0;
97721e72a00SBarry Smith   while (ak < an && bk < bn) {
97821e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
97921e72a00SBarry Smith       L_[k] = aI[ak];
98021e72a00SBarry Smith       ++ak;
98121e72a00SBarry Smith       ++bk;
98221e72a00SBarry Smith       ++k;
98321e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
98421e72a00SBarry Smith       L_[k] = aI[ak];
98521e72a00SBarry Smith       ++ak;
98621e72a00SBarry Smith       ++k;
98721e72a00SBarry Smith     } else {
98821e72a00SBarry Smith       L_[k] = bI[bk];
98921e72a00SBarry Smith       ++bk;
99021e72a00SBarry Smith       ++k;
99121e72a00SBarry Smith     }
99221e72a00SBarry Smith   }
99321e72a00SBarry Smith   if (ak < an) {
9949566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
99521e72a00SBarry Smith     k += (an - ak);
99621e72a00SBarry Smith   }
99721e72a00SBarry Smith   if (bk < bn) {
9989566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
99921e72a00SBarry Smith     k += (bn - bk);
100021e72a00SBarry Smith   }
100121e72a00SBarry Smith   *n = k;
10023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100321e72a00SBarry Smith }
1004b4301105SBarry Smith 
1005c1f0200aSDmitry Karpeev /*@
1006811af0c4SBarry Smith   PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
1007c1f0200aSDmitry Karpeev   The additional arrays are the same length as sorted arrays and are merged
1008c1f0200aSDmitry Karpeev   in the order determined by the merging of the sorted pair.
1009c1f0200aSDmitry Karpeev 
1010c1f0200aSDmitry Karpeev   Not Collective
1011c1f0200aSDmitry Karpeev 
1012c1f0200aSDmitry Karpeev   Input Parameters:
1013c1f0200aSDmitry Karpeev + an - number of values in the first array
10143334a463SJames Wright . aI - first sorted array of `PetscInt`
10153334a463SJames Wright . aJ - first additional array of `PetscInt`
1016c1f0200aSDmitry Karpeev . bn - number of values in the second array
10173334a463SJames Wright . bI - second array of `PetscInt`
10183334a463SJames Wright - bJ - second additional of `PetscInt`
1019c1f0200aSDmitry Karpeev 
1020c1f0200aSDmitry Karpeev   Output Parameters:
1021c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn)
102214c49006SJed Brown . L - merged sorted array
1023c1f0200aSDmitry Karpeev - J - merged additional array
1024c1f0200aSDmitry Karpeev 
1025811af0c4SBarry Smith   Note:
1026da81f932SPierre 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
1027811af0c4SBarry Smith 
1028c1f0200aSDmitry Karpeev   Level: intermediate
1029c1f0200aSDmitry Karpeev 
1030db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1031c1f0200aSDmitry Karpeev @*/
1032d71ae5a4SJacob 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)
1033d71ae5a4SJacob Faibussowitsch {
103470d8d27cSBarry Smith   PetscInt n_, *L_, *J_, ak, bk, k;
1035c1f0200aSDmitry Karpeev 
103614c49006SJed Brown   PetscFunctionBegin;
10374f572ea9SToby Isaac   PetscAssertPointer(L, 8);
10384f572ea9SToby Isaac   PetscAssertPointer(J, 9);
1039c1f0200aSDmitry Karpeev   n_ = an + bn;
1040c1f0200aSDmitry Karpeev   *n = n_;
104148a46eb9SPierre Jolivet   if (!*L) PetscCall(PetscMalloc1(n_, L));
1042d7aa01a8SHong Zhang   L_ = *L;
104348a46eb9SPierre Jolivet   if (!*J) PetscCall(PetscMalloc1(n_, J));
1044c1f0200aSDmitry Karpeev   J_ = *J;
1045c1f0200aSDmitry Karpeev   k = ak = bk = 0;
1046c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
1047c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
1048d7aa01a8SHong Zhang       L_[k] = aI[ak];
1049c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
1050c1f0200aSDmitry Karpeev       ++ak;
1051c1f0200aSDmitry Karpeev       ++k;
1052a297a907SKarl Rupp     } else {
1053d7aa01a8SHong Zhang       L_[k] = bI[bk];
1054c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
1055c1f0200aSDmitry Karpeev       ++bk;
1056c1f0200aSDmitry Karpeev       ++k;
1057c1f0200aSDmitry Karpeev     }
1058c1f0200aSDmitry Karpeev   }
1059c1f0200aSDmitry Karpeev   if (ak < an) {
10609566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
10619566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1062c1f0200aSDmitry Karpeev     k += (an - ak);
1063c1f0200aSDmitry Karpeev   }
1064c1f0200aSDmitry Karpeev   if (bk < bn) {
10659566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
10669566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1067c1f0200aSDmitry Karpeev   }
10683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1069c1f0200aSDmitry Karpeev }
1070c1f0200aSDmitry Karpeev 
1071e498c390SJed Brown /*@
1072811af0c4SBarry Smith   PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
1073e498c390SJed Brown 
1074e498c390SJed Brown   Not Collective
1075e498c390SJed Brown 
1076e498c390SJed Brown   Input Parameters:
1077e498c390SJed Brown + an - number of values in the first array
10783334a463SJames Wright . aI - first sorted array of `PetscMPIInt`
1079e498c390SJed Brown . bn - number of values in the second array
10803334a463SJames Wright - bI - second array of `PetscMPIInt`
1081e498c390SJed Brown 
1082e498c390SJed Brown   Output Parameters:
1083e498c390SJed Brown + n - number of values in the merged array (<= an + bn)
1084e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed
1085e498c390SJed Brown 
1086e498c390SJed Brown   Level: intermediate
1087e498c390SJed Brown 
1088db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1089e498c390SJed Brown @*/
1090d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1091d71ae5a4SJacob Faibussowitsch {
1092e498c390SJed Brown   PetscInt ai, bi, k;
1093e498c390SJed Brown 
1094e498c390SJed Brown   PetscFunctionBegin;
109557508eceSPierre Jolivet   if (!*L) PetscCall(PetscMalloc1(an + bn, L));
1096e498c390SJed Brown   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1097e498c390SJed Brown     PetscInt t = -1;
1098e498c390SJed Brown     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
1099fbccb6d4SPierre Jolivet     for (; bi < bn && bI[bi] == t; bi++);
1100e498c390SJed Brown     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
1101fbccb6d4SPierre Jolivet     for (; ai < an && aI[ai] == t; ai++);
1102e498c390SJed Brown   }
1103e498c390SJed Brown   *n = k;
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1105e498c390SJed Brown }
1106e5c89e4eSSatish Balay 
11076c2863d0SBarry Smith /*@C
110842eaa7f3SBarry Smith   PetscProcessTree - Prepares tree data to be displayed graphically
110942eaa7f3SBarry Smith 
1110cc4c1da9SBarry Smith   Not Collective, No Fortran Support
111142eaa7f3SBarry Smith 
111242eaa7f3SBarry Smith   Input Parameters:
111342eaa7f3SBarry Smith + n        - number of values
111442eaa7f3SBarry Smith . mask     - indicates those entries in the tree, location 0 is always masked
111542eaa7f3SBarry Smith - parentid - indicates the parent of each entry
111642eaa7f3SBarry Smith 
111742eaa7f3SBarry Smith   Output Parameters:
111842eaa7f3SBarry Smith + Nlevels   - the number of levels
111942eaa7f3SBarry Smith . Level     - for each node tells its level
1120aec76313SJacob Faibussowitsch . Levelcnt  - the number of nodes on each level
112142eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc
112242eaa7f3SBarry Smith - Column    - for each id tells its column index
112342eaa7f3SBarry Smith 
11246c2863d0SBarry Smith   Level: developer
112542eaa7f3SBarry Smith 
1126811af0c4SBarry Smith   Note:
112795452b02SPatrick Sanan   This code is not currently used
112842eaa7f3SBarry Smith 
1129db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
113042eaa7f3SBarry Smith @*/
1131d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1132d71ae5a4SJacob Faibussowitsch {
113342eaa7f3SBarry Smith   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1134ace3abfcSBarry Smith   PetscBool done = PETSC_FALSE;
113542eaa7f3SBarry Smith 
113642eaa7f3SBarry Smith   PetscFunctionBegin;
113708401ef6SPierre Jolivet   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
113842eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
113942eaa7f3SBarry Smith     if (mask[i]) continue;
114008401ef6SPierre Jolivet     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1141cc73adaaSBarry Smith     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
114242eaa7f3SBarry Smith   }
114342eaa7f3SBarry Smith 
114442eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
114542eaa7f3SBarry Smith     if (!mask[i]) nmask++;
114642eaa7f3SBarry Smith   }
114742eaa7f3SBarry Smith 
114842eaa7f3SBarry Smith   /* determine the level in the tree of each node */
11499566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &level));
1150a297a907SKarl Rupp 
115142eaa7f3SBarry Smith   level[0] = 1;
115242eaa7f3SBarry Smith   while (!done) {
115342eaa7f3SBarry Smith     done = PETSC_TRUE;
115442eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
115542eaa7f3SBarry Smith       if (mask[i]) continue;
115642eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
115742eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
115842eaa7f3SBarry Smith     }
115942eaa7f3SBarry Smith   }
116042eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
116142eaa7f3SBarry Smith     level[i]--;
116242eaa7f3SBarry Smith     nlevels = PetscMax(nlevels, level[i]);
116342eaa7f3SBarry Smith   }
116442eaa7f3SBarry Smith 
116542eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
11669566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nlevels, &levelcnt));
116742eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
116842eaa7f3SBarry Smith     if (mask[i]) continue;
116942eaa7f3SBarry Smith     levelcnt[level[i] - 1]++;
117042eaa7f3SBarry Smith   }
1171a297a907SKarl Rupp   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
117242eaa7f3SBarry Smith 
117342eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
11749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
11759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmask, &idbylevel));
117642eaa7f3SBarry Smith   for (j = 1; j <= nlevels; j++) {
117742eaa7f3SBarry Smith     cnt = 0;
117842eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
117942eaa7f3SBarry Smith       if (mask[i]) continue;
118042eaa7f3SBarry Smith       if (level[i] != j) continue;
118142eaa7f3SBarry Smith       workid[cnt]         = i;
118242eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
118342eaa7f3SBarry Smith     }
118442eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
118542eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
11869566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
118742eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
118842eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
11899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
119042eaa7f3SBarry Smith     tcnt += cnt;
119142eaa7f3SBarry Smith   }
119208401ef6SPierre Jolivet   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
11939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(workid, workparentid));
119442eaa7f3SBarry Smith 
119542eaa7f3SBarry Smith   /* for each node list its column */
11969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &column));
119742eaa7f3SBarry Smith   cnt = 0;
119842eaa7f3SBarry Smith   for (j = 0; j < nlevels; j++) {
1199ad540459SPierre Jolivet     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
120042eaa7f3SBarry Smith   }
120142eaa7f3SBarry Smith 
120242eaa7f3SBarry Smith   *Nlevels   = nlevels;
120342eaa7f3SBarry Smith   *Level     = level;
120442eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
120542eaa7f3SBarry Smith   *Idbylevel = idbylevel;
120642eaa7f3SBarry Smith   *Column    = column;
12073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120842eaa7f3SBarry Smith }
1209ce605777SToby Isaac 
1210ce605777SToby Isaac /*@
1211811af0c4SBarry Smith   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1212ce605777SToby Isaac 
1213ce605777SToby Isaac   Collective
1214ce605777SToby Isaac 
1215ce605777SToby Isaac   Input Parameters:
1216ce605777SToby Isaac + comm - the MPI communicator
12173334a463SJames Wright . n    - the local number of `PetscInt`
12183334a463SJames Wright - keys - the local array of `PetscInt`
1219ce605777SToby Isaac 
1220ce605777SToby Isaac   Output Parameters:
1221ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1222ce605777SToby Isaac 
1223ce605777SToby Isaac   Level: developer
1224ce605777SToby Isaac 
1225db781477SPatrick Sanan .seealso: `PetscParallelSortInt()`
1226ce605777SToby Isaac @*/
1227d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1228d71ae5a4SJacob Faibussowitsch {
1229ce605777SToby Isaac   PetscBool   sorted;
1230ce605777SToby Isaac   PetscInt    i, min, max, prevmax;
12312a1da528SToby Isaac   PetscMPIInt rank;
1232ce605777SToby Isaac 
1233ce605777SToby Isaac   PetscFunctionBegin;
1234ce605777SToby Isaac   sorted = PETSC_TRUE;
12351690c2aeSBarry Smith   min    = PETSC_INT_MAX;
12361690c2aeSBarry Smith   max    = PETSC_INT_MIN;
1237ce605777SToby Isaac   if (n) {
1238ce605777SToby Isaac     min = keys[0];
1239ce605777SToby Isaac     max = keys[0];
1240ce605777SToby Isaac   }
1241ce605777SToby Isaac   for (i = 1; i < n; i++) {
1242ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1243ce605777SToby Isaac     min = PetscMin(min, keys[i]);
1244ce605777SToby Isaac     max = PetscMax(max, keys[i]);
1245ce605777SToby Isaac   }
1246ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
12471690c2aeSBarry Smith   prevmax = PETSC_INT_MIN;
12489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
12499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
12501690c2aeSBarry Smith   if (rank == 0) prevmax = PETSC_INT_MIN;
1251ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1252462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
12533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1254ce605777SToby Isaac }
1255