xref: /petsc/src/sys/utils/sorti.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
4c4762a1bSJed Brown    One can use src/sys/tests/ex52.c for benchmarking.
5e5c89e4eSSatish Balay  */
6af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"  I*/
7f1cab4e1SJunchao Zhang #include <petsc/private/hashseti.h>
8e5c89e4eSSatish Balay 
99371c9d4SSatish Balay #define MEDIAN3(v, a, b, c) (v[a] < v[b] ? (v[b] < v[c] ? (b) : (v[a] < v[c] ? (c) : (a))) : (v[c] < v[b] ? (b) : (v[a] < v[c] ? (a) : (c))))
10a8582168SJed Brown 
11a8582168SJed Brown #define MEDIAN(v, right) MEDIAN3(v, right / 4, right / 2, right / 4 * 3)
12ef8e3583SJed Brown 
1359e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */
149371c9d4SSatish Balay #define SWAP1(a, b, t1) \
159371c9d4SSatish Balay   do { \
169371c9d4SSatish Balay     t1 = a; \
179371c9d4SSatish Balay     a  = b; \
189371c9d4SSatish Balay     b  = t1; \
199371c9d4SSatish Balay   } while (0)
209371c9d4SSatish Balay #define SWAP2(a, b, c, d, t1, t2) \
219371c9d4SSatish Balay   do { \
229371c9d4SSatish Balay     t1 = a; \
239371c9d4SSatish Balay     a  = b; \
249371c9d4SSatish Balay     b  = t1; \
259371c9d4SSatish Balay     t2 = c; \
269371c9d4SSatish Balay     c  = d; \
279371c9d4SSatish Balay     d  = t2; \
289371c9d4SSatish Balay   } while (0)
299371c9d4SSatish Balay #define SWAP3(a, b, c, d, e, f, t1, t2, t3) \
309371c9d4SSatish Balay   do { \
319371c9d4SSatish Balay     t1 = a; \
329371c9d4SSatish Balay     a  = b; \
339371c9d4SSatish Balay     b  = t1; \
349371c9d4SSatish Balay     t2 = c; \
359371c9d4SSatish Balay     c  = d; \
369371c9d4SSatish Balay     d  = t2; \
379371c9d4SSatish Balay     t3 = e; \
389371c9d4SSatish Balay     e  = f; \
399371c9d4SSatish Balay     f  = t3; \
409371c9d4SSatish Balay   } while (0)
4159e16d97SJunchao Zhang 
4259e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
4359e16d97SJunchao Zhang #define SWAP2Data(a, b, c, d, t1, t2, siz) \
4459e16d97SJunchao Zhang   do { \
459371c9d4SSatish Balay     t1 = a; \
469371c9d4SSatish Balay     a  = b; \
479371c9d4SSatish Balay     b  = t1; \
489566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(t2, c, siz)); \
499566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(c, d, siz)); \
509566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(d, t2, siz)); \
5159e16d97SJunchao Zhang   } while (0)
52e5c89e4eSSatish Balay 
53e5c89e4eSSatish Balay /*
5459e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
55e5c89e4eSSatish Balay 
5659e16d97SJunchao Zhang    Input Parameters:
5759e16d97SJunchao Zhang     + X         - array to partition
5859e16d97SJunchao Zhang     . pivot     - a pivot of X[]
5959e16d97SJunchao Zhang     . t1        - temp variable for X
6059e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
6159e16d97SJunchao Zhang 
6259e16d97SJunchao Zhang    Output Parameters:
6359e16d97SJunchao Zhang     + l,r       - of type PetscInt
6459e16d97SJunchao Zhang 
65811af0c4SBarry Smith    Note:
6659e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
6759e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
6859e16d97SJunchao Zhang  */
6959e16d97SJunchao Zhang #define TwoWayPartition1(X, pivot, t1, lo, hi, l, r) \
7059e16d97SJunchao Zhang   do { \
7159e16d97SJunchao Zhang     l = lo; \
7259e16d97SJunchao Zhang     r = hi; \
7359e16d97SJunchao Zhang     while (1) { \
7459e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
7559e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
769371c9d4SSatish Balay       if (l >= r) { \
779371c9d4SSatish Balay         r++; \
789371c9d4SSatish Balay         break; \
799371c9d4SSatish Balay       } \
8059e16d97SJunchao Zhang       SWAP1(X[l], X[r], t1); \
8159e16d97SJunchao Zhang       l++; \
8259e16d97SJunchao Zhang       r--; \
8359e16d97SJunchao Zhang     } \
8459e16d97SJunchao Zhang   } while (0)
8559e16d97SJunchao Zhang 
86ce605777SToby Isaac /*
87ce605777SToby Isaac    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
88ce605777SToby Isaac 
89ce605777SToby Isaac    Input Parameters:
90ce605777SToby Isaac     + X         - array to partition
91ce605777SToby Isaac     . pivot     - a pivot of X[]
92ce605777SToby Isaac     . t1        - temp variable for X
93ce605777SToby Isaac     - lo,hi     - lower and upper bound of the array
94ce605777SToby Isaac 
95ce605777SToby Isaac    Output Parameters:
96ce605777SToby Isaac     + l,r       - of type PetscInt
97ce605777SToby Isaac 
98811af0c4SBarry Smith    Note:
99ce605777SToby Isaac     The TwoWayPartition2/3 variants also partition other arrays along with X.
100ce605777SToby Isaac     These arrays can have different types, so they provide their own temp t2,t3
101ce605777SToby Isaac  */
102ce605777SToby Isaac #define TwoWayPartitionReverse1(X, pivot, t1, lo, hi, l, r) \
103ce605777SToby Isaac   do { \
104ce605777SToby Isaac     l = lo; \
105ce605777SToby Isaac     r = hi; \
106ce605777SToby Isaac     while (1) { \
107ce605777SToby Isaac       while (X[l] > pivot) l++; \
108ce605777SToby Isaac       while (X[r] < pivot) r--; \
1099371c9d4SSatish Balay       if (l >= r) { \
1109371c9d4SSatish Balay         r++; \
1119371c9d4SSatish Balay         break; \
1129371c9d4SSatish Balay       } \
113ce605777SToby Isaac       SWAP1(X[l], X[r], t1); \
114ce605777SToby Isaac       l++; \
115ce605777SToby Isaac       r--; \
116ce605777SToby Isaac     } \
117ce605777SToby Isaac   } while (0)
118ce605777SToby Isaac 
11959e16d97SJunchao Zhang #define TwoWayPartition2(X, Y, pivot, t1, t2, lo, hi, l, r) \
12059e16d97SJunchao Zhang   do { \
12159e16d97SJunchao Zhang     l = lo; \
12259e16d97SJunchao Zhang     r = hi; \
12359e16d97SJunchao Zhang     while (1) { \
12459e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
12559e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1269371c9d4SSatish Balay       if (l >= r) { \
1279371c9d4SSatish Balay         r++; \
1289371c9d4SSatish Balay         break; \
1299371c9d4SSatish Balay       } \
13059e16d97SJunchao Zhang       SWAP2(X[l], X[r], Y[l], Y[r], t1, t2); \
13159e16d97SJunchao Zhang       l++; \
13259e16d97SJunchao Zhang       r--; \
13359e16d97SJunchao Zhang     } \
13459e16d97SJunchao Zhang   } while (0)
13559e16d97SJunchao Zhang 
13659e16d97SJunchao Zhang #define TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, lo, hi, l, r) \
13759e16d97SJunchao Zhang   do { \
13859e16d97SJunchao Zhang     l = lo; \
13959e16d97SJunchao Zhang     r = hi; \
14059e16d97SJunchao Zhang     while (1) { \
14159e16d97SJunchao Zhang       while (X[l] < pivot) l++; \
14259e16d97SJunchao Zhang       while (X[r] > pivot) r--; \
1439371c9d4SSatish Balay       if (l >= r) { \
1449371c9d4SSatish Balay         r++; \
1459371c9d4SSatish Balay         break; \
1469371c9d4SSatish Balay       } \
14759e16d97SJunchao Zhang       SWAP3(X[l], X[r], Y[l], Y[r], Z[l], Z[r], t1, t2, t3); \
14859e16d97SJunchao Zhang       l++; \
14959e16d97SJunchao Zhang       r--; \
15059e16d97SJunchao Zhang     } \
15159e16d97SJunchao Zhang   } while (0)
15259e16d97SJunchao Zhang 
15359e16d97SJunchao Zhang /* Templates for similar functions used below */
1545f80ce2aSJacob Faibussowitsch #define QuickSort1(FuncName, X, n, pivot, t1) \
15559e16d97SJunchao Zhang   do { \
156981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
15759e16d97SJunchao Zhang     if (n < 8) { \
15859e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
15959e16d97SJunchao Zhang         pivot = X[i]; \
16059e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
16159e16d97SJunchao Zhang           if (pivot > X[j]) { \
16259e16d97SJunchao Zhang             SWAP1(X[i], X[j], t1); \
16359e16d97SJunchao Zhang             pivot = X[i]; \
16459e16d97SJunchao Zhang           } \
16559e16d97SJunchao Zhang         } \
16659e16d97SJunchao Zhang       } \
16759e16d97SJunchao Zhang     } else { \
16859e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
16959e16d97SJunchao Zhang       pivot = X[p]; \
17059e16d97SJunchao Zhang       TwoWayPartition1(X, pivot, t1, 0, hi, l, r); \
1719566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1729566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
17359e16d97SJunchao Zhang     } \
17459e16d97SJunchao Zhang   } while (0)
17559e16d97SJunchao Zhang 
176ce605777SToby Isaac /* Templates for similar functions used below */
1775f80ce2aSJacob Faibussowitsch #define QuickSortReverse1(FuncName, X, n, pivot, t1) \
178ce605777SToby Isaac   do { \
179981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
180ce605777SToby Isaac     if (n < 8) { \
181ce605777SToby Isaac       for (i = 0; i < n; i++) { \
182ce605777SToby Isaac         pivot = X[i]; \
183ce605777SToby Isaac         for (j = i + 1; j < n; j++) { \
184ce605777SToby Isaac           if (pivot < X[j]) { \
185ce605777SToby Isaac             SWAP1(X[i], X[j], t1); \
186ce605777SToby Isaac             pivot = X[i]; \
187ce605777SToby Isaac           } \
188ce605777SToby Isaac         } \
189ce605777SToby Isaac       } \
190ce605777SToby Isaac     } else { \
191ce605777SToby Isaac       p     = MEDIAN(X, hi); \
192ce605777SToby Isaac       pivot = X[p]; \
193ce605777SToby Isaac       TwoWayPartitionReverse1(X, pivot, t1, 0, hi, l, r); \
1949566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X)); \
1959566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r)); \
196ce605777SToby Isaac     } \
197ce605777SToby Isaac   } while (0)
198ce605777SToby Isaac 
1995f80ce2aSJacob Faibussowitsch #define QuickSort2(FuncName, X, Y, n, pivot, t1, t2) \
20059e16d97SJunchao Zhang   do { \
201981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
20259e16d97SJunchao Zhang     if (n < 8) { \
20359e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
20459e16d97SJunchao Zhang         pivot = X[i]; \
20559e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
20659e16d97SJunchao Zhang           if (pivot > X[j]) { \
20759e16d97SJunchao Zhang             SWAP2(X[i], X[j], Y[i], Y[j], t1, t2); \
20859e16d97SJunchao Zhang             pivot = X[i]; \
20959e16d97SJunchao Zhang           } \
21059e16d97SJunchao Zhang         } \
21159e16d97SJunchao Zhang       } \
21259e16d97SJunchao Zhang     } else { \
21359e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
21459e16d97SJunchao Zhang       pivot = X[p]; \
21559e16d97SJunchao Zhang       TwoWayPartition2(X, Y, pivot, t1, t2, 0, hi, l, r); \
2169566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y)); \
2179566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r)); \
21859e16d97SJunchao Zhang     } \
21959e16d97SJunchao Zhang   } while (0)
22059e16d97SJunchao Zhang 
2215f80ce2aSJacob Faibussowitsch #define QuickSort3(FuncName, X, Y, Z, n, pivot, t1, t2, t3) \
22259e16d97SJunchao Zhang   do { \
223981bb840SJunchao Zhang     PetscCount i, j, p, l, r, hi = n - 1; \
22459e16d97SJunchao Zhang     if (n < 8) { \
22559e16d97SJunchao Zhang       for (i = 0; i < n; i++) { \
22659e16d97SJunchao Zhang         pivot = X[i]; \
22759e16d97SJunchao Zhang         for (j = i + 1; j < n; j++) { \
22859e16d97SJunchao Zhang           if (pivot > X[j]) { \
22959e16d97SJunchao Zhang             SWAP3(X[i], X[j], Y[i], Y[j], Z[i], Z[j], t1, t2, t3); \
23059e16d97SJunchao Zhang             pivot = X[i]; \
23159e16d97SJunchao Zhang           } \
23259e16d97SJunchao Zhang         } \
23359e16d97SJunchao Zhang       } \
23459e16d97SJunchao Zhang     } else { \
23559e16d97SJunchao Zhang       p     = MEDIAN(X, hi); \
23659e16d97SJunchao Zhang       pivot = X[p]; \
23759e16d97SJunchao Zhang       TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, 0, hi, l, r); \
2389566063dSJacob Faibussowitsch       PetscCall(FuncName(l, X, Y, Z)); \
2399566063dSJacob Faibussowitsch       PetscCall(FuncName(hi - r + 1, X + r, Y + r, Z + r)); \
24059e16d97SJunchao Zhang     } \
24159e16d97SJunchao Zhang   } while (0)
242e5c89e4eSSatish Balay 
243e5c89e4eSSatish Balay /*@
244811af0c4SBarry Smith    PetscSortedInt - Determines whether the `PetscInt` array is sorted.
2456a7c706bSVaclav Hapla 
2466a7c706bSVaclav Hapla    Not Collective
2476a7c706bSVaclav Hapla 
2486a7c706bSVaclav Hapla    Input Parameters:
2496a7c706bSVaclav Hapla +  n  - number of values
2506a7c706bSVaclav Hapla -  X  - array of integers
2516a7c706bSVaclav Hapla 
2526a7c706bSVaclav Hapla    Output Parameters:
2536a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
2546a7c706bSVaclav Hapla 
2556a7c706bSVaclav Hapla    Level: intermediate
2566a7c706bSVaclav Hapla 
257db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
2586a7c706bSVaclav Hapla @*/
259*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedInt(PetscInt n, const PetscInt X[], PetscBool *sorted)
260*d71ae5a4SJacob Faibussowitsch {
2616a7c706bSVaclav Hapla   PetscFunctionBegin;
2625f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
2635f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(sorted, 3);
2646a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
2656a7c706bSVaclav Hapla   PetscFunctionReturn(0);
2666a7c706bSVaclav Hapla }
2676a7c706bSVaclav Hapla 
2686a7c706bSVaclav Hapla /*@
269811af0c4SBarry Smith    PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.
270e5c89e4eSSatish Balay 
271e5c89e4eSSatish Balay    Not Collective
272e5c89e4eSSatish Balay 
273e5c89e4eSSatish Balay    Input Parameters:
274e5c89e4eSSatish Balay +  n  - number of values
27559e16d97SJunchao Zhang -  X  - array of integers
276e5c89e4eSSatish Balay 
277811af0c4SBarry Smith    Note:
278811af0c4SBarry Smith    This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
279a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
280676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
281676f2a66SJacob Faibussowitsch 
282e5c89e4eSSatish Balay    Level: intermediate
283e5c89e4eSSatish Balay 
284db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
285e5c89e4eSSatish Balay @*/
286*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[])
287*d71ae5a4SJacob Faibussowitsch {
28859e16d97SJunchao Zhang   PetscInt pivot, t1;
289e5c89e4eSSatish Balay 
290e5c89e4eSSatish Balay   PetscFunctionBegin;
2915f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
2925f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortInt, X, n, pivot, t1);
293e5c89e4eSSatish Balay   PetscFunctionReturn(0);
294e5c89e4eSSatish Balay }
295e5c89e4eSSatish Balay 
2961db36a52SBarry Smith /*@
297811af0c4SBarry Smith    PetscSortReverseInt - Sorts an array of `PetscInt` in place in decreasing order.
298ce605777SToby Isaac 
299ce605777SToby Isaac    Not Collective
300ce605777SToby Isaac 
301ce605777SToby Isaac    Input Parameters:
302ce605777SToby Isaac +  n  - number of values
303ce605777SToby Isaac -  X  - array of integers
304ce605777SToby Isaac 
305ce605777SToby Isaac    Level: intermediate
306ce605777SToby Isaac 
307db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
308ce605777SToby Isaac @*/
309*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[])
310*d71ae5a4SJacob Faibussowitsch {
311ce605777SToby Isaac   PetscInt pivot, t1;
312ce605777SToby Isaac 
313ce605777SToby Isaac   PetscFunctionBegin;
3145f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
3155f80ce2aSJacob Faibussowitsch   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
316ce605777SToby Isaac   PetscFunctionReturn(0);
317ce605777SToby Isaac }
318ce605777SToby Isaac 
319ce605777SToby Isaac /*@
320811af0c4SBarry Smith    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array
32122ab5688SLisandro Dalcin 
32222ab5688SLisandro Dalcin    Not Collective
32322ab5688SLisandro Dalcin 
32422ab5688SLisandro Dalcin    Input Parameters:
32522ab5688SLisandro Dalcin +  n  - number of values
32659e16d97SJunchao Zhang -  X  - sorted array of integers
32722ab5688SLisandro Dalcin 
32822ab5688SLisandro Dalcin    Output Parameter:
32922ab5688SLisandro Dalcin .  n - number of non-redundant values
33022ab5688SLisandro Dalcin 
33122ab5688SLisandro Dalcin    Level: intermediate
33222ab5688SLisandro Dalcin 
333db781477SPatrick Sanan .seealso: `PetscSortInt()`
33422ab5688SLisandro Dalcin @*/
335*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[])
336*d71ae5a4SJacob Faibussowitsch {
33722ab5688SLisandro Dalcin   PetscInt i, s = 0, N = *n, b = 0;
33822ab5688SLisandro Dalcin 
33922ab5688SLisandro Dalcin   PetscFunctionBegin;
3405f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n, 1);
341fb61b9e4SStefano Zampini   PetscCheckSorted(*n, X);
34222ab5688SLisandro Dalcin   for (i = 0; i < N - 1; i++) {
34359e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
3449371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
3459371c9d4SSatish Balay       b++;
34622ab5688SLisandro Dalcin     } else s++;
34722ab5688SLisandro Dalcin   }
34822ab5688SLisandro Dalcin   *n = N - s;
34922ab5688SLisandro Dalcin   PetscFunctionReturn(0);
35022ab5688SLisandro Dalcin }
35122ab5688SLisandro Dalcin 
35222ab5688SLisandro Dalcin /*@
353811af0c4SBarry Smith    PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates
354b6a18d21SVaclav Hapla 
355b6a18d21SVaclav Hapla    Not Collective
356b6a18d21SVaclav Hapla 
357b6a18d21SVaclav Hapla    Input Parameters:
358b6a18d21SVaclav Hapla +  n  - number of values
359b6a18d21SVaclav Hapla -  X  - sorted array of integers
360b6a18d21SVaclav Hapla 
361b6a18d21SVaclav Hapla    Output Parameter:
362b6a18d21SVaclav Hapla .  dups - True if the array has dups, otherwise false
363b6a18d21SVaclav Hapla 
364b6a18d21SVaclav Hapla    Level: intermediate
365b6a18d21SVaclav Hapla 
366db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
367b6a18d21SVaclav Hapla @*/
368*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg)
369*d71ae5a4SJacob Faibussowitsch {
370b6a18d21SVaclav Hapla   PetscInt i;
371b6a18d21SVaclav Hapla 
372b6a18d21SVaclav Hapla   PetscFunctionBegin;
373b6a18d21SVaclav Hapla   PetscCheckSorted(n, X);
374b6a18d21SVaclav Hapla   *flg = PETSC_FALSE;
375b6a18d21SVaclav Hapla   for (i = 0; i < n - 1; i++) {
376b6a18d21SVaclav Hapla     if (X[i + 1] == X[i]) {
377b6a18d21SVaclav Hapla       *flg = PETSC_TRUE;
378b6a18d21SVaclav Hapla       break;
379b6a18d21SVaclav Hapla     }
380b6a18d21SVaclav Hapla   }
381b6a18d21SVaclav Hapla   PetscFunctionReturn(0);
382b6a18d21SVaclav Hapla }
383b6a18d21SVaclav Hapla 
384b6a18d21SVaclav Hapla /*@
385811af0c4SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
3861db36a52SBarry Smith 
3871db36a52SBarry Smith    Not Collective
3881db36a52SBarry Smith 
3891db36a52SBarry Smith    Input Parameters:
3901db36a52SBarry Smith +  n  - number of values
39159e16d97SJunchao Zhang -  X  - array of integers
3921db36a52SBarry Smith 
3931db36a52SBarry Smith    Output Parameter:
3941db36a52SBarry Smith .  n - number of non-redundant values
3951db36a52SBarry Smith 
3961db36a52SBarry Smith    Level: intermediate
3971db36a52SBarry Smith 
398db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
3991db36a52SBarry Smith @*/
400*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
401*d71ae5a4SJacob Faibussowitsch {
4021db36a52SBarry Smith   PetscFunctionBegin;
4035f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n, 1);
4049566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(*n, X));
4059566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(n, X));
4061db36a52SBarry Smith   PetscFunctionReturn(0);
4071db36a52SBarry Smith }
4081db36a52SBarry Smith 
40960e03357SMatthew G Knepley /*@
410811af0c4SBarry Smith  PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt`
41160e03357SMatthew G Knepley 
41260e03357SMatthew G Knepley    Not Collective
41360e03357SMatthew G Knepley 
41460e03357SMatthew G Knepley    Input Parameters:
41560e03357SMatthew G Knepley +  key - the integer to locate
416d9f1162dSJed Brown .  n   - number of values in the array
41759e16d97SJunchao Zhang -  X  - array of integers
41860e03357SMatthew G Knepley 
41960e03357SMatthew G Knepley    Output Parameter:
420d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
42160e03357SMatthew G Knepley 
42260e03357SMatthew G Knepley    Level: intermediate
42360e03357SMatthew G Knepley 
424db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
42560e03357SMatthew G Knepley @*/
426*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
427*d71ae5a4SJacob Faibussowitsch {
428d9f1162dSJed Brown   PetscInt lo = 0, hi = n;
42960e03357SMatthew G Knepley 
43060e03357SMatthew G Knepley   PetscFunctionBegin;
4315f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(loc, 4);
4329371c9d4SSatish Balay   if (!n) {
4339371c9d4SSatish Balay     *loc = -1;
4349371c9d4SSatish Balay     PetscFunctionReturn(0);
4359371c9d4SSatish Balay   }
436dadcf809SJacob Faibussowitsch   PetscValidIntPointer(X, 3);
4376a7c706bSVaclav Hapla   PetscCheckSorted(n, X);
438d9f1162dSJed Brown   while (hi - lo > 1) {
439d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo) / 2;
44059e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
441d9f1162dSJed Brown     else lo = mid;
44260e03357SMatthew G Knepley   }
44359e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
44460e03357SMatthew G Knepley   PetscFunctionReturn(0);
44560e03357SMatthew G Knepley }
44660e03357SMatthew G Knepley 
447d2aeb606SJed Brown /*@
448811af0c4SBarry Smith   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
449f1cab4e1SJunchao Zhang 
450f1cab4e1SJunchao Zhang    Not Collective
451f1cab4e1SJunchao Zhang 
452f1cab4e1SJunchao Zhang    Input Parameters:
453f1cab4e1SJunchao Zhang +  n  - number of values in the array
454f1cab4e1SJunchao Zhang -  X  - array of integers
455f1cab4e1SJunchao Zhang 
456f1cab4e1SJunchao Zhang    Output Parameter:
457f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
458f1cab4e1SJunchao Zhang 
459f1cab4e1SJunchao Zhang    Level: intermediate
460f1cab4e1SJunchao Zhang 
461db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
462f1cab4e1SJunchao Zhang @*/
463*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
464*d71ae5a4SJacob Faibussowitsch {
465f1cab4e1SJunchao Zhang   PetscInt   i;
466f1cab4e1SJunchao Zhang   PetscHSetI ht;
467f1cab4e1SJunchao Zhang   PetscBool  missing;
468f1cab4e1SJunchao Zhang 
469f1cab4e1SJunchao Zhang   PetscFunctionBegin;
4705f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X, 2);
4715f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(dups, 3);
472f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
473f1cab4e1SJunchao Zhang   if (n > 1) {
4749566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&ht));
4759566063dSJacob Faibussowitsch     PetscCall(PetscHSetIResize(ht, n));
476f1cab4e1SJunchao Zhang     for (i = 0; i < n; i++) {
4779566063dSJacob Faibussowitsch       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
4789371c9d4SSatish Balay       if (!missing) {
4799371c9d4SSatish Balay         *dups = PETSC_TRUE;
4809371c9d4SSatish Balay         break;
4819371c9d4SSatish Balay       }
482f1cab4e1SJunchao Zhang     }
4839566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&ht));
484f1cab4e1SJunchao Zhang   }
485f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
486f1cab4e1SJunchao Zhang }
487f1cab4e1SJunchao Zhang 
488f1cab4e1SJunchao Zhang /*@
489811af0c4SBarry Smith   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
490d2aeb606SJed Brown 
491d2aeb606SJed Brown    Not Collective
492d2aeb606SJed Brown 
493d2aeb606SJed Brown    Input Parameters:
494d2aeb606SJed Brown +  key - the integer to locate
495d2aeb606SJed Brown .  n   - number of values in the array
49659e16d97SJunchao Zhang -  X   - array of integers
497d2aeb606SJed Brown 
498d2aeb606SJed Brown    Output Parameter:
499d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
500d2aeb606SJed Brown 
501d2aeb606SJed Brown    Level: intermediate
502d2aeb606SJed Brown 
503db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
504d2aeb606SJed Brown @*/
505*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
506*d71ae5a4SJacob Faibussowitsch {
507d2aeb606SJed Brown   PetscInt lo = 0, hi = n;
508d2aeb606SJed Brown 
509d2aeb606SJed Brown   PetscFunctionBegin;
510dadcf809SJacob Faibussowitsch   PetscValidIntPointer(loc, 4);
5119371c9d4SSatish Balay   if (!n) {
5129371c9d4SSatish Balay     *loc = -1;
5139371c9d4SSatish Balay     PetscFunctionReturn(0);
5149371c9d4SSatish Balay   }
515dadcf809SJacob Faibussowitsch   PetscValidIntPointer(X, 3);
5166a7c706bSVaclav Hapla   PetscCheckSorted(n, X);
517d2aeb606SJed Brown   while (hi - lo > 1) {
518d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo) / 2;
51959e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
520d2aeb606SJed Brown     else lo = mid;
521d2aeb606SJed Brown   }
52259e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
523e5c89e4eSSatish Balay   PetscFunctionReturn(0);
524e5c89e4eSSatish Balay }
525e5c89e4eSSatish Balay 
526e5c89e4eSSatish Balay /*@
527811af0c4SBarry Smith    PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
528811af0c4SBarry Smith        changes a second array of `PetscInt` to match the sorted first array.
529e5c89e4eSSatish Balay 
530e5c89e4eSSatish Balay    Not Collective
531e5c89e4eSSatish Balay 
532e5c89e4eSSatish Balay    Input Parameters:
533e5c89e4eSSatish Balay +  n  - number of values
53459e16d97SJunchao Zhang .  X  - array of integers
53559e16d97SJunchao Zhang -  Y  - second array of integers
536e5c89e4eSSatish Balay 
537e5c89e4eSSatish Balay    Level: intermediate
538e5c89e4eSSatish Balay 
539db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
540e5c89e4eSSatish Balay @*/
541*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[])
542*d71ae5a4SJacob Faibussowitsch {
54359e16d97SJunchao Zhang   PetscInt pivot, t1, t2;
544e5c89e4eSSatish Balay 
545e5c89e4eSSatish Balay   PetscFunctionBegin;
5465f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
547c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
548c1f0200aSDmitry Karpeev }
549c1f0200aSDmitry Karpeev 
550c1f0200aSDmitry Karpeev /*@
551811af0c4SBarry Smith    PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
552811af0c4SBarry Smith        changes a pair of `PetscInt` arrays to match the sorted first array.
553c1f0200aSDmitry Karpeev 
554c1f0200aSDmitry Karpeev    Not Collective
555c1f0200aSDmitry Karpeev 
556c1f0200aSDmitry Karpeev    Input Parameters:
557c1f0200aSDmitry Karpeev +  n  - number of values
55859e16d97SJunchao Zhang .  X  - array of integers
55959e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
56059e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
561c1f0200aSDmitry Karpeev 
562c1f0200aSDmitry Karpeev    Level: intermediate
563c1f0200aSDmitry Karpeev 
564db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
565c1f0200aSDmitry Karpeev @*/
566*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[])
567*d71ae5a4SJacob Faibussowitsch {
56859e16d97SJunchao Zhang   PetscInt pivot, t1, t2, t3;
569c1f0200aSDmitry Karpeev 
570c1f0200aSDmitry Karpeev   PetscFunctionBegin;
5715f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
572c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
573c1f0200aSDmitry Karpeev }
574c1f0200aSDmitry Karpeev 
57517d7d925SStefano Zampini /*@
576811af0c4SBarry Smith    PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
577811af0c4SBarry Smith        changes a second array of `PetscCount` to match the sorted first array.
578981bb840SJunchao Zhang 
579981bb840SJunchao Zhang    Not Collective
580981bb840SJunchao Zhang 
581981bb840SJunchao Zhang    Input Parameters:
582981bb840SJunchao Zhang +  n  - number of values
583981bb840SJunchao Zhang .  X  - array of integers
584981bb840SJunchao Zhang -  Y  - second array of PetscCounts (signed integers)
585981bb840SJunchao Zhang 
586981bb840SJunchao Zhang    Level: intermediate
587981bb840SJunchao Zhang 
588db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
589981bb840SJunchao Zhang @*/
590*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
591*d71ae5a4SJacob Faibussowitsch {
592981bb840SJunchao Zhang   PetscInt   pivot, t1;
593981bb840SJunchao Zhang   PetscCount t2;
594981bb840SJunchao Zhang 
595981bb840SJunchao Zhang   PetscFunctionBegin;
5965f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
597981bb840SJunchao Zhang   PetscFunctionReturn(0);
598981bb840SJunchao Zhang }
599981bb840SJunchao Zhang 
600981bb840SJunchao Zhang /*@
601811af0c4SBarry Smith    PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
602811af0c4SBarry Smith        changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.
603981bb840SJunchao Zhang 
604981bb840SJunchao Zhang    Not Collective
605981bb840SJunchao Zhang 
606981bb840SJunchao Zhang    Input Parameters:
607981bb840SJunchao Zhang +  n  - number of values
608981bb840SJunchao Zhang .  X  - array of integers
609981bb840SJunchao Zhang .  Y  - second array of integers (first array of the pair)
610981bb840SJunchao Zhang -  Z  - third array of PetscCounts  (second array of the pair)
611981bb840SJunchao Zhang 
612981bb840SJunchao Zhang    Level: intermediate
613981bb840SJunchao Zhang 
614811af0c4SBarry Smith    Note:
615981bb840SJunchao 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.
616981bb840SJunchao Zhang 
617db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
618981bb840SJunchao Zhang @*/
619*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
620*d71ae5a4SJacob Faibussowitsch {
621981bb840SJunchao Zhang   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
622981bb840SJunchao Zhang   PetscCount t3;            /* temp for Z[] */
623981bb840SJunchao Zhang 
624981bb840SJunchao Zhang   PetscFunctionBegin;
6255f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
626981bb840SJunchao Zhang   PetscFunctionReturn(0);
627981bb840SJunchao Zhang }
628981bb840SJunchao Zhang 
629981bb840SJunchao Zhang /*@
630811af0c4SBarry Smith   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
6316a7c706bSVaclav Hapla 
6326a7c706bSVaclav Hapla    Not Collective
6336a7c706bSVaclav Hapla 
6346a7c706bSVaclav Hapla    Input Parameters:
6356a7c706bSVaclav Hapla +  n  - number of values
6366a7c706bSVaclav Hapla -  X  - array of integers
6376a7c706bSVaclav Hapla 
6386a7c706bSVaclav Hapla    Output Parameters:
6396a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
6406a7c706bSVaclav Hapla 
6416a7c706bSVaclav Hapla    Level: intermediate
6426a7c706bSVaclav Hapla 
643db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
6446a7c706bSVaclav Hapla @*/
645*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted)
646*d71ae5a4SJacob Faibussowitsch {
6476a7c706bSVaclav Hapla   PetscFunctionBegin;
6486a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
6496a7c706bSVaclav Hapla   PetscFunctionReturn(0);
6506a7c706bSVaclav Hapla }
6516a7c706bSVaclav Hapla 
6526a7c706bSVaclav Hapla /*@
653811af0c4SBarry Smith    PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
65417d7d925SStefano Zampini 
65517d7d925SStefano Zampini    Not Collective
65617d7d925SStefano Zampini 
65717d7d925SStefano Zampini    Input Parameters:
65817d7d925SStefano Zampini +  n  - number of values
65959e16d97SJunchao Zhang -  X  - array of integers
66017d7d925SStefano Zampini 
66117d7d925SStefano Zampini    Level: intermediate
66217d7d925SStefano Zampini 
663811af0c4SBarry Smith    Note:
664676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
665a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
666676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
667676f2a66SJacob Faibussowitsch 
668db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
66917d7d925SStefano Zampini @*/
670*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[])
671*d71ae5a4SJacob Faibussowitsch {
67259e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
67317d7d925SStefano Zampini 
67417d7d925SStefano Zampini   PetscFunctionBegin;
6755f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
67617d7d925SStefano Zampini   PetscFunctionReturn(0);
67717d7d925SStefano Zampini }
67817d7d925SStefano Zampini 
67917d7d925SStefano Zampini /*@
680811af0c4SBarry Smith    PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
68117d7d925SStefano Zampini 
68217d7d925SStefano Zampini    Not Collective
68317d7d925SStefano Zampini 
68417d7d925SStefano Zampini    Input Parameters:
68517d7d925SStefano Zampini +  n  - number of values
68659e16d97SJunchao Zhang -  X  - array of integers
68717d7d925SStefano Zampini 
68817d7d925SStefano Zampini    Output Parameter:
68917d7d925SStefano Zampini .  n - number of non-redundant values
69017d7d925SStefano Zampini 
69117d7d925SStefano Zampini    Level: intermediate
69217d7d925SStefano Zampini 
693db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
69417d7d925SStefano Zampini @*/
695*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
696*d71ae5a4SJacob Faibussowitsch {
6975f80ce2aSJacob Faibussowitsch   PetscInt s = 0, N = *n, b = 0;
69817d7d925SStefano Zampini 
69917d7d925SStefano Zampini   PetscFunctionBegin;
7009566063dSJacob Faibussowitsch   PetscCall(PetscSortMPIInt(N, X));
7015f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < N - 1; i++) {
70259e16d97SJunchao Zhang     if (X[b + s + 1] != X[b]) {
7039371c9d4SSatish Balay       X[b + 1] = X[b + s + 1];
7049371c9d4SSatish Balay       b++;
705a297a907SKarl Rupp     } else s++;
70617d7d925SStefano Zampini   }
70717d7d925SStefano Zampini   *n = N - s;
70817d7d925SStefano Zampini   PetscFunctionReturn(0);
70917d7d925SStefano Zampini }
710c1f0200aSDmitry Karpeev 
7114d615ea0SBarry Smith /*@
712811af0c4SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
713811af0c4SBarry Smith        changes a second `PetscMPIInt` array to match the sorted first array.
7144d615ea0SBarry Smith 
7154d615ea0SBarry Smith    Not Collective
7164d615ea0SBarry Smith 
7174d615ea0SBarry Smith    Input Parameters:
7184d615ea0SBarry Smith +  n  - number of values
71959e16d97SJunchao Zhang .  X  - array of integers
72059e16d97SJunchao Zhang -  Y  - second array of integers
7214d615ea0SBarry Smith 
7224d615ea0SBarry Smith    Level: intermediate
7234d615ea0SBarry Smith 
724db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
7254d615ea0SBarry Smith @*/
726*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[])
727*d71ae5a4SJacob Faibussowitsch {
72859e16d97SJunchao Zhang   PetscMPIInt pivot, t1, t2;
7294d615ea0SBarry Smith 
7304d615ea0SBarry Smith   PetscFunctionBegin;
7315f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
7325569a785SJunchao Zhang   PetscFunctionReturn(0);
7335569a785SJunchao Zhang }
7345569a785SJunchao Zhang 
7355569a785SJunchao Zhang /*@
736811af0c4SBarry Smith    PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
737811af0c4SBarry Smith        changes a second array of `PetscInt` to match the sorted first array.
7385569a785SJunchao Zhang 
7395569a785SJunchao Zhang    Not Collective
7405569a785SJunchao Zhang 
7415569a785SJunchao Zhang    Input Parameters:
7425569a785SJunchao Zhang +  n  - number of values
74359e16d97SJunchao Zhang .  X  - array of MPI integers
74459e16d97SJunchao Zhang -  Y  - second array of Petsc integers
7455569a785SJunchao Zhang 
7465569a785SJunchao Zhang    Level: intermediate
7475569a785SJunchao Zhang 
748811af0c4SBarry Smith    Note:
749811af0c4SBarry Smith    This routine is useful when one needs to sort MPI ranks with other integer arrays.
7505569a785SJunchao Zhang 
751db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
7525569a785SJunchao Zhang @*/
753*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[])
754*d71ae5a4SJacob Faibussowitsch {
75559e16d97SJunchao Zhang   PetscMPIInt pivot, t1;
7565569a785SJunchao Zhang   PetscInt    t2;
7575569a785SJunchao Zhang 
7585569a785SJunchao Zhang   PetscFunctionBegin;
7595f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
760e5c89e4eSSatish Balay   PetscFunctionReturn(0);
761e5c89e4eSSatish Balay }
762e5c89e4eSSatish Balay 
763e5c89e4eSSatish Balay /*@
764811af0c4SBarry Smith    PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
765811af0c4SBarry Smith        changes a second `PetscScalar` array to match the sorted first array.
766e5c89e4eSSatish Balay 
767e5c89e4eSSatish Balay    Not Collective
768e5c89e4eSSatish Balay 
769e5c89e4eSSatish Balay    Input Parameters:
770e5c89e4eSSatish Balay +  n  - number of values
77159e16d97SJunchao Zhang .  X  - array of integers
77259e16d97SJunchao Zhang -  Y  - second array of scalars
773e5c89e4eSSatish Balay 
774e5c89e4eSSatish Balay    Level: intermediate
775e5c89e4eSSatish Balay 
776db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
777e5c89e4eSSatish Balay @*/
778*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[])
779*d71ae5a4SJacob Faibussowitsch {
78059e16d97SJunchao Zhang   PetscInt    pivot, t1;
78159e16d97SJunchao Zhang   PetscScalar t2;
782e5c89e4eSSatish Balay 
783e5c89e4eSSatish Balay   PetscFunctionBegin;
7845f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
78517df18f2SToby Isaac   PetscFunctionReturn(0);
78617df18f2SToby Isaac }
78717df18f2SToby Isaac 
7886c2863d0SBarry Smith /*@C
789811af0c4SBarry Smith    PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
79017df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
79117df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
79217df18f2SToby Isaac 
79317df18f2SToby Isaac    Not Collective
79417df18f2SToby Isaac 
79517df18f2SToby Isaac    Input Parameters:
79617df18f2SToby Isaac +  n  - number of values
79759e16d97SJunchao Zhang .  X  - array of integers
79859e16d97SJunchao Zhang .  Y  - second array of data
79917df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
80059e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
80117df18f2SToby Isaac 
80217df18f2SToby Isaac    Level: intermediate
80317df18f2SToby Isaac 
804db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
80517df18f2SToby Isaac @*/
806*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2)
807*d71ae5a4SJacob Faibussowitsch {
80859e16d97SJunchao Zhang   char    *YY = (char *)Y;
8095f80ce2aSJacob Faibussowitsch   PetscInt t1, pivot, hi = n - 1;
81017df18f2SToby Isaac 
81117df18f2SToby Isaac   PetscFunctionBegin;
81217df18f2SToby Isaac   if (n < 8) {
8135f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < n; i++) {
81459e16d97SJunchao Zhang       pivot = X[i];
8155f80ce2aSJacob Faibussowitsch       for (PetscInt j = i + 1; j < n; j++) {
81659e16d97SJunchao Zhang         if (pivot > X[j]) {
81759e16d97SJunchao Zhang           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
81859e16d97SJunchao Zhang           pivot = X[i];
81917df18f2SToby Isaac         }
82017df18f2SToby Isaac       }
82117df18f2SToby Isaac     }
82217df18f2SToby Isaac   } else {
82359e16d97SJunchao Zhang     /* Two way partition */
8245f80ce2aSJacob Faibussowitsch     PetscInt l = 0, r = hi;
8255f80ce2aSJacob Faibussowitsch 
8265f80ce2aSJacob Faibussowitsch     pivot = X[MEDIAN(X, hi)];
82759e16d97SJunchao Zhang     while (1) {
82859e16d97SJunchao Zhang       while (X[l] < pivot) l++;
82959e16d97SJunchao Zhang       while (X[r] > pivot) r--;
8309371c9d4SSatish Balay       if (l >= r) {
8319371c9d4SSatish Balay         r++;
8329371c9d4SSatish Balay         break;
8339371c9d4SSatish Balay       }
83459e16d97SJunchao Zhang       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
83559e16d97SJunchao Zhang       l++;
83659e16d97SJunchao Zhang       r--;
83759e16d97SJunchao Zhang     }
8389566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
8399566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
84017df18f2SToby Isaac   }
84117df18f2SToby Isaac   PetscFunctionReturn(0);
84217df18f2SToby Isaac }
84317df18f2SToby Isaac 
84421e72a00SBarry Smith /*@
845811af0c4SBarry Smith    PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.
84621e72a00SBarry Smith 
84721e72a00SBarry Smith    Not Collective
84821e72a00SBarry Smith 
84921e72a00SBarry Smith    Input Parameters:
85021e72a00SBarry Smith +  an  - number of values in the first array
85121e72a00SBarry Smith .  aI  - first sorted array of integers
85221e72a00SBarry Smith .  bn  - number of values in the second array
85321e72a00SBarry Smith -  bI  - second array of integers
85421e72a00SBarry Smith 
85521e72a00SBarry Smith    Output Parameters:
85621e72a00SBarry Smith +  n   - number of values in the merged array
8576c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
85821e72a00SBarry Smith 
85921e72a00SBarry Smith    Level: intermediate
86021e72a00SBarry Smith 
861db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
86221e72a00SBarry Smith @*/
863*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
864*d71ae5a4SJacob Faibussowitsch {
86521e72a00SBarry Smith   PetscInt *L_ = *L, ak, bk, k;
86621e72a00SBarry Smith 
867362febeeSStefano Zampini   PetscFunctionBegin;
86821e72a00SBarry Smith   if (!L_) {
8699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(an + bn, L));
87021e72a00SBarry Smith     L_ = *L;
87121e72a00SBarry Smith   }
87221e72a00SBarry Smith   k = ak = bk = 0;
87321e72a00SBarry Smith   while (ak < an && bk < bn) {
87421e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
87521e72a00SBarry Smith       L_[k] = aI[ak];
87621e72a00SBarry Smith       ++ak;
87721e72a00SBarry Smith       ++bk;
87821e72a00SBarry Smith       ++k;
87921e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
88021e72a00SBarry Smith       L_[k] = aI[ak];
88121e72a00SBarry Smith       ++ak;
88221e72a00SBarry Smith       ++k;
88321e72a00SBarry Smith     } else {
88421e72a00SBarry Smith       L_[k] = bI[bk];
88521e72a00SBarry Smith       ++bk;
88621e72a00SBarry Smith       ++k;
88721e72a00SBarry Smith     }
88821e72a00SBarry Smith   }
88921e72a00SBarry Smith   if (ak < an) {
8909566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
89121e72a00SBarry Smith     k += (an - ak);
89221e72a00SBarry Smith   }
89321e72a00SBarry Smith   if (bk < bn) {
8949566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
89521e72a00SBarry Smith     k += (bn - bk);
89621e72a00SBarry Smith   }
89721e72a00SBarry Smith   *n = k;
89821e72a00SBarry Smith   PetscFunctionReturn(0);
89921e72a00SBarry Smith }
900b4301105SBarry Smith 
901c1f0200aSDmitry Karpeev /*@
902811af0c4SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
903c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
904c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
905c1f0200aSDmitry Karpeev 
906c1f0200aSDmitry Karpeev    Not Collective
907c1f0200aSDmitry Karpeev 
908c1f0200aSDmitry Karpeev    Input Parameters:
909c1f0200aSDmitry Karpeev +  an  - number of values in the first array
910c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
911c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
912c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
913c1f0200aSDmitry Karpeev .  bI  - second array of integers
914c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
915c1f0200aSDmitry Karpeev 
916c1f0200aSDmitry Karpeev    Output Parameters:
917c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
91814c49006SJed Brown .  L   - merged sorted array
919c1f0200aSDmitry Karpeev -  J   - merged additional array
920c1f0200aSDmitry Karpeev 
921811af0c4SBarry Smith    Note:
92295452b02SPatrick Sanan     if L or J point to non-null arrays then this routine will assume they are of the approproate size and use them, otherwise this routine will allocate space for them
923811af0c4SBarry Smith 
924c1f0200aSDmitry Karpeev    Level: intermediate
925c1f0200aSDmitry Karpeev 
926db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
927c1f0200aSDmitry Karpeev @*/
928*d71ae5a4SJacob 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)
929*d71ae5a4SJacob Faibussowitsch {
93070d8d27cSBarry Smith   PetscInt n_, *L_, *J_, ak, bk, k;
931c1f0200aSDmitry Karpeev 
93214c49006SJed Brown   PetscFunctionBegin;
933dadcf809SJacob Faibussowitsch   PetscValidPointer(L, 8);
934dadcf809SJacob Faibussowitsch   PetscValidPointer(J, 9);
935c1f0200aSDmitry Karpeev   n_ = an + bn;
936c1f0200aSDmitry Karpeev   *n = n_;
93748a46eb9SPierre Jolivet   if (!*L) PetscCall(PetscMalloc1(n_, L));
938d7aa01a8SHong Zhang   L_ = *L;
93948a46eb9SPierre Jolivet   if (!*J) PetscCall(PetscMalloc1(n_, J));
940c1f0200aSDmitry Karpeev   J_ = *J;
941c1f0200aSDmitry Karpeev   k = ak = bk = 0;
942c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
943c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
944d7aa01a8SHong Zhang       L_[k] = aI[ak];
945c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
946c1f0200aSDmitry Karpeev       ++ak;
947c1f0200aSDmitry Karpeev       ++k;
948a297a907SKarl Rupp     } else {
949d7aa01a8SHong Zhang       L_[k] = bI[bk];
950c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
951c1f0200aSDmitry Karpeev       ++bk;
952c1f0200aSDmitry Karpeev       ++k;
953c1f0200aSDmitry Karpeev     }
954c1f0200aSDmitry Karpeev   }
955c1f0200aSDmitry Karpeev   if (ak < an) {
9569566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
9579566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
958c1f0200aSDmitry Karpeev     k += (an - ak);
959c1f0200aSDmitry Karpeev   }
960c1f0200aSDmitry Karpeev   if (bk < bn) {
9619566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
9629566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
963c1f0200aSDmitry Karpeev   }
964c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
965c1f0200aSDmitry Karpeev }
966c1f0200aSDmitry Karpeev 
967e498c390SJed Brown /*@
968811af0c4SBarry Smith    PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.
969e498c390SJed Brown 
970e498c390SJed Brown    Not Collective
971e498c390SJed Brown 
972e498c390SJed Brown    Input Parameters:
973e498c390SJed Brown +  an  - number of values in the first array
974e498c390SJed Brown .  aI  - first sorted array of integers
975e498c390SJed Brown .  bn  - number of values in the second array
976e498c390SJed Brown -  bI  - second array of integers
977e498c390SJed Brown 
978e498c390SJed Brown    Output Parameters:
979e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
980e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
981e498c390SJed Brown 
982e498c390SJed Brown    Level: intermediate
983e498c390SJed Brown 
984db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
985e498c390SJed Brown @*/
986*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
987*d71ae5a4SJacob Faibussowitsch {
988e498c390SJed Brown   PetscInt ai, bi, k;
989e498c390SJed Brown 
990e498c390SJed Brown   PetscFunctionBegin;
9919566063dSJacob Faibussowitsch   if (!*L) PetscCall(PetscMalloc1((an + bn), L));
992e498c390SJed Brown   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
993e498c390SJed Brown     PetscInt t = -1;
994e498c390SJed Brown     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
9959371c9d4SSatish Balay     for (; bi < bn && bI[bi] == t; bi++)
9969371c9d4SSatish Balay       ;
997e498c390SJed Brown     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
9989371c9d4SSatish Balay     for (; ai < an && aI[ai] == t; ai++)
9999371c9d4SSatish Balay       ;
1000e498c390SJed Brown   }
1001e498c390SJed Brown   *n = k;
1002e498c390SJed Brown   PetscFunctionReturn(0);
1003e498c390SJed Brown }
1004e5c89e4eSSatish Balay 
10056c2863d0SBarry Smith /*@C
100642eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
100742eaa7f3SBarry Smith 
100842eaa7f3SBarry Smith    Not Collective
100942eaa7f3SBarry Smith 
101042eaa7f3SBarry Smith    Input Parameters:
101142eaa7f3SBarry Smith +  n  - number of values
101242eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
101342eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
101442eaa7f3SBarry Smith 
101542eaa7f3SBarry Smith    Output Parameters:
101642eaa7f3SBarry Smith +  Nlevels - the number of levels
101742eaa7f3SBarry Smith .  Level - for each node tells its level
101842eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
101942eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
102042eaa7f3SBarry Smith -  Column - for each id tells its column index
102142eaa7f3SBarry Smith 
10226c2863d0SBarry Smith    Level: developer
102342eaa7f3SBarry Smith 
1024811af0c4SBarry Smith    Note:
102595452b02SPatrick Sanan     This code is not currently used
102642eaa7f3SBarry Smith 
1027db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
102842eaa7f3SBarry Smith @*/
1029*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1030*d71ae5a4SJacob Faibussowitsch {
103142eaa7f3SBarry Smith   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1032ace3abfcSBarry Smith   PetscBool done = PETSC_FALSE;
103342eaa7f3SBarry Smith 
103442eaa7f3SBarry Smith   PetscFunctionBegin;
103508401ef6SPierre Jolivet   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
103642eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
103742eaa7f3SBarry Smith     if (mask[i]) continue;
103808401ef6SPierre Jolivet     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1039cc73adaaSBarry Smith     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
104042eaa7f3SBarry Smith   }
104142eaa7f3SBarry Smith 
104242eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
104342eaa7f3SBarry Smith     if (!mask[i]) nmask++;
104442eaa7f3SBarry Smith   }
104542eaa7f3SBarry Smith 
104642eaa7f3SBarry Smith   /* determine the level in the tree of each node */
10479566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &level));
1048a297a907SKarl Rupp 
104942eaa7f3SBarry Smith   level[0] = 1;
105042eaa7f3SBarry Smith   while (!done) {
105142eaa7f3SBarry Smith     done = PETSC_TRUE;
105242eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
105342eaa7f3SBarry Smith       if (mask[i]) continue;
105442eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
105542eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
105642eaa7f3SBarry Smith     }
105742eaa7f3SBarry Smith   }
105842eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
105942eaa7f3SBarry Smith     level[i]--;
106042eaa7f3SBarry Smith     nlevels = PetscMax(nlevels, level[i]);
106142eaa7f3SBarry Smith   }
106242eaa7f3SBarry Smith 
106342eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
10649566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nlevels, &levelcnt));
106542eaa7f3SBarry Smith   for (i = 0; i < n; i++) {
106642eaa7f3SBarry Smith     if (mask[i]) continue;
106742eaa7f3SBarry Smith     levelcnt[level[i] - 1]++;
106842eaa7f3SBarry Smith   }
1069a297a907SKarl Rupp   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
107042eaa7f3SBarry Smith 
107142eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
10729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
10739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmask, &idbylevel));
107442eaa7f3SBarry Smith   for (j = 1; j <= nlevels; j++) {
107542eaa7f3SBarry Smith     cnt = 0;
107642eaa7f3SBarry Smith     for (i = 0; i < n; i++) {
107742eaa7f3SBarry Smith       if (mask[i]) continue;
107842eaa7f3SBarry Smith       if (level[i] != j) continue;
107942eaa7f3SBarry Smith       workid[cnt]         = i;
108042eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
108142eaa7f3SBarry Smith     }
108242eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
108342eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
10849566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
108542eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
108642eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
10879566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
108842eaa7f3SBarry Smith     tcnt += cnt;
108942eaa7f3SBarry Smith   }
109008401ef6SPierre Jolivet   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
10919566063dSJacob Faibussowitsch   PetscCall(PetscFree2(workid, workparentid));
109242eaa7f3SBarry Smith 
109342eaa7f3SBarry Smith   /* for each node list its column */
10949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &column));
109542eaa7f3SBarry Smith   cnt = 0;
109642eaa7f3SBarry Smith   for (j = 0; j < nlevels; j++) {
1097ad540459SPierre Jolivet     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
109842eaa7f3SBarry Smith   }
109942eaa7f3SBarry Smith 
110042eaa7f3SBarry Smith   *Nlevels   = nlevels;
110142eaa7f3SBarry Smith   *Level     = level;
110242eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
110342eaa7f3SBarry Smith   *Idbylevel = idbylevel;
110442eaa7f3SBarry Smith   *Column    = column;
110542eaa7f3SBarry Smith   PetscFunctionReturn(0);
110642eaa7f3SBarry Smith }
1107ce605777SToby Isaac 
1108ce605777SToby Isaac /*@
1109811af0c4SBarry Smith   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1110ce605777SToby Isaac 
1111ce605777SToby Isaac   Collective
1112ce605777SToby Isaac 
1113ce605777SToby Isaac   Input Parameters:
1114ce605777SToby Isaac + comm - the MPI communicator
1115ce605777SToby Isaac . n - the local number of integers
1116ce605777SToby Isaac - keys - the local array of integers
1117ce605777SToby Isaac 
1118ce605777SToby Isaac   Output Parameters:
1119ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1120ce605777SToby Isaac 
1121ce605777SToby Isaac   Level: developer
1122ce605777SToby Isaac 
1123db781477SPatrick Sanan .seealso: `PetscParallelSortInt()`
1124ce605777SToby Isaac @*/
1125*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1126*d71ae5a4SJacob Faibussowitsch {
1127ce605777SToby Isaac   PetscBool   sorted;
1128ce605777SToby Isaac   PetscInt    i, min, max, prevmax;
11292a1da528SToby Isaac   PetscMPIInt rank;
1130ce605777SToby Isaac 
1131ce605777SToby Isaac   PetscFunctionBegin;
1132ce605777SToby Isaac   sorted = PETSC_TRUE;
1133ce605777SToby Isaac   min    = PETSC_MAX_INT;
1134ce605777SToby Isaac   max    = PETSC_MIN_INT;
1135ce605777SToby Isaac   if (n) {
1136ce605777SToby Isaac     min = keys[0];
1137ce605777SToby Isaac     max = keys[0];
1138ce605777SToby Isaac   }
1139ce605777SToby Isaac   for (i = 1; i < n; i++) {
1140ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1141ce605777SToby Isaac     min = PetscMin(min, keys[i]);
1142ce605777SToby Isaac     max = PetscMax(max, keys[i]);
1143ce605777SToby Isaac   }
1144ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1145ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
11469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
11479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1148dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1149ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
11509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
1151ce605777SToby Isaac   PetscFunctionReturn(0);
1152ce605777SToby Isaac }
1153