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 @*/
PetscSortedInt(PetscCount n,const PetscInt X[],PetscBool * sorted)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 @*/
PetscSortedInt64(PetscCount n,const PetscInt64 X[],PetscBool * sorted)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 @*/
PetscSortInt(PetscCount n,PetscInt X[])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 @*/
PetscSortInt64(PetscCount n,PetscInt64 X[])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 @*/
PetscSortCount(PetscCount n,PetscCount X[])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 @*/
PetscSortReverseInt(PetscCount n,PetscInt X[])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 @*/
PetscSortedRemoveDupsInt(PetscInt * n,PetscInt X[])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 @*/
PetscSortedCheckDupsInt(PetscCount n,const PetscInt X[],PetscBool * flg)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 /*@
46277aa4f85SJames Wright PetscSortedCheckDupsCount - Checks if a sorted `PetscCount` array has duplicates
46377aa4f85SJames Wright
46477aa4f85SJames Wright Not Collective
46577aa4f85SJames Wright
46677aa4f85SJames Wright Input Parameters:
46777aa4f85SJames Wright + n - number of values
46877aa4f85SJames Wright - X - sorted array of `PetscCount`
46977aa4f85SJames Wright
47077aa4f85SJames Wright Output Parameter:
47177aa4f85SJames Wright . flg - True if the array has duplications, otherwise false
47277aa4f85SJames Wright
47377aa4f85SJames Wright Level: intermediate
47477aa4f85SJames Wright
47598480730SJames Wright .seealso: `PetscCount`, `PetscSortCount()`, `PetscSortedCheckDupsInt()`
47677aa4f85SJames Wright @*/
PetscSortedCheckDupsCount(PetscCount n,const PetscCount X[],PetscBool * flg)47777aa4f85SJames Wright PetscErrorCode PetscSortedCheckDupsCount(PetscCount n, const PetscCount X[], PetscBool *flg)
47877aa4f85SJames Wright {
47977aa4f85SJames Wright PetscCount i;
48077aa4f85SJames Wright
48177aa4f85SJames Wright PetscFunctionBegin;
48277aa4f85SJames Wright PetscCheckSorted(n, X);
48377aa4f85SJames Wright *flg = PETSC_FALSE;
48477aa4f85SJames Wright for (i = 0; i < n - 1; i++) {
48577aa4f85SJames Wright if (X[i + 1] == X[i]) {
48677aa4f85SJames Wright *flg = PETSC_TRUE;
48777aa4f85SJames Wright break;
48877aa4f85SJames Wright }
48977aa4f85SJames Wright }
49077aa4f85SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
49177aa4f85SJames Wright }
49277aa4f85SJames Wright
49377aa4f85SJames Wright /*@
494811af0c4SBarry Smith PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries
4951db36a52SBarry Smith
4961db36a52SBarry Smith Not Collective
4971db36a52SBarry Smith
4981db36a52SBarry Smith Input Parameters:
4991db36a52SBarry Smith + n - number of values
5003334a463SJames Wright - X - array of `PetscInt`
5011db36a52SBarry Smith
5021db36a52SBarry Smith Output Parameter:
5031db36a52SBarry Smith . n - number of non-redundant values
5041db36a52SBarry Smith
5051db36a52SBarry Smith Level: intermediate
5061db36a52SBarry Smith
507db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
5081db36a52SBarry Smith @*/
PetscSortRemoveDupsInt(PetscInt * n,PetscInt X[])509d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
510d71ae5a4SJacob Faibussowitsch {
5111db36a52SBarry Smith PetscFunctionBegin;
5124f572ea9SToby Isaac PetscAssertPointer(n, 1);
5139566063dSJacob Faibussowitsch PetscCall(PetscSortInt(*n, X));
5149566063dSJacob Faibussowitsch PetscCall(PetscSortedRemoveDupsInt(n, X));
5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5161db36a52SBarry Smith }
5171db36a52SBarry Smith
51860e03357SMatthew G Knepley /*@
5193334a463SJames Wright PetscFindInt - Finds the location of a `PetscInt` key in a sorted array of `PetscInt`
52060e03357SMatthew G Knepley
52160e03357SMatthew G Knepley Not Collective
52260e03357SMatthew G Knepley
52360e03357SMatthew G Knepley Input Parameters:
5243334a463SJames Wright + key - the `PetscInt` key to locate
525d9f1162dSJed Brown . n - number of values in the array
5263334a463SJames Wright - X - array of `PetscInt`
52760e03357SMatthew G Knepley
52860e03357SMatthew G Knepley Output Parameter:
529d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
53060e03357SMatthew G Knepley
53160e03357SMatthew G Knepley Level: intermediate
53260e03357SMatthew G Knepley
533db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
53460e03357SMatthew G Knepley @*/
PetscFindInt(PetscInt key,PetscCount n,const PetscInt X[],PetscInt * loc)5356497c311SBarry Smith PetscErrorCode PetscFindInt(PetscInt key, PetscCount n, const PetscInt X[], PetscInt *loc)
536d71ae5a4SJacob Faibussowitsch {
5376497c311SBarry Smith PetscInt lo = 0, hi;
53860e03357SMatthew G Knepley
53960e03357SMatthew G Knepley PetscFunctionBegin;
5404f572ea9SToby Isaac PetscAssertPointer(loc, 4);
5419371c9d4SSatish Balay if (!n) {
5429371c9d4SSatish Balay *loc = -1;
5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5449371c9d4SSatish Balay }
5454f572ea9SToby Isaac PetscAssertPointer(X, 3);
5466497c311SBarry Smith PetscCall(PetscIntCast(n, &hi));
547d9f1162dSJed Brown while (hi - lo > 1) {
548d9f1162dSJed Brown PetscInt mid = lo + (hi - lo) / 2;
54942f965a0SMatthew 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]);
55059e16d97SJunchao Zhang if (key < X[mid]) hi = mid;
551d9f1162dSJed Brown else lo = mid;
55260e03357SMatthew G Knepley }
55359e16d97SJunchao Zhang *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55560e03357SMatthew G Knepley }
55660e03357SMatthew G Knepley
557d2aeb606SJed Brown /*@
55877aa4f85SJames Wright PetscFindCount - Finds the location of a `PetscCount` key in a sorted array of `PetscCount`
55977aa4f85SJames Wright
56077aa4f85SJames Wright Not Collective
56177aa4f85SJames Wright
56277aa4f85SJames Wright Input Parameters:
56377aa4f85SJames Wright + key - the `PetscCount` key to locate
56477aa4f85SJames Wright . n - number of values in the array
56577aa4f85SJames Wright - X - array of `PetscCount`
56677aa4f85SJames Wright
56777aa4f85SJames Wright Output Parameter:
56877aa4f85SJames Wright . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
56977aa4f85SJames Wright
57077aa4f85SJames Wright Level: intermediate
57177aa4f85SJames Wright
57298480730SJames Wright .seealso: `PetscCount`, `PetscSortCount()`
57377aa4f85SJames Wright @*/
PetscFindCount(PetscCount key,PetscCount n,const PetscCount X[],PetscCount * loc)57477aa4f85SJames Wright PetscErrorCode PetscFindCount(PetscCount key, PetscCount n, const PetscCount X[], PetscCount *loc)
57577aa4f85SJames Wright {
57677aa4f85SJames Wright PetscCount lo = 0, hi;
57777aa4f85SJames Wright
57877aa4f85SJames Wright PetscFunctionBegin;
57977aa4f85SJames Wright PetscAssertPointer(loc, 4);
58077aa4f85SJames Wright if (!n) {
58177aa4f85SJames Wright *loc = -1;
58277aa4f85SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
58377aa4f85SJames Wright }
58477aa4f85SJames Wright PetscAssertPointer(X, 3);
58577aa4f85SJames Wright hi = n;
58677aa4f85SJames Wright while (hi - lo > 1) {
58777aa4f85SJames Wright PetscCount mid = lo + (hi - lo) / 2;
58877aa4f85SJames Wright PetscAssert(X[lo] <= X[mid] && X[mid] <= X[hi - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array was not sorted: (%" PetscCount_FMT ", %" PetscCount_FMT ", %" PetscCount_FMT ")", X[lo], X[mid], X[hi - 1]);
58977aa4f85SJames Wright if (key < X[mid]) hi = mid;
59077aa4f85SJames Wright else lo = mid;
59177aa4f85SJames Wright }
59277aa4f85SJames Wright *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
59377aa4f85SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
59477aa4f85SJames Wright }
59577aa4f85SJames Wright
59677aa4f85SJames Wright /*@
597811af0c4SBarry Smith PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates
598f1cab4e1SJunchao Zhang
599f1cab4e1SJunchao Zhang Not Collective
600f1cab4e1SJunchao Zhang
601f1cab4e1SJunchao Zhang Input Parameters:
602f1cab4e1SJunchao Zhang + n - number of values in the array
6033334a463SJames Wright - X - array of `PetscInt`
604f1cab4e1SJunchao Zhang
605f1cab4e1SJunchao Zhang Output Parameter:
606f1cab4e1SJunchao Zhang . dups - True if the array has dups, otherwise false
607f1cab4e1SJunchao Zhang
608f1cab4e1SJunchao Zhang Level: intermediate
609f1cab4e1SJunchao Zhang
610db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
611f1cab4e1SJunchao Zhang @*/
PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool * dups)612d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
613d71ae5a4SJacob Faibussowitsch {
614f1cab4e1SJunchao Zhang PetscInt i;
615f1cab4e1SJunchao Zhang PetscHSetI ht;
616f1cab4e1SJunchao Zhang PetscBool missing;
617f1cab4e1SJunchao Zhang
618f1cab4e1SJunchao Zhang PetscFunctionBegin;
6194f572ea9SToby Isaac if (n) PetscAssertPointer(X, 2);
6204f572ea9SToby Isaac PetscAssertPointer(dups, 3);
621f1cab4e1SJunchao Zhang *dups = PETSC_FALSE;
622f1cab4e1SJunchao Zhang if (n > 1) {
6239566063dSJacob Faibussowitsch PetscCall(PetscHSetICreate(&ht));
6249566063dSJacob Faibussowitsch PetscCall(PetscHSetIResize(ht, n));
625f1cab4e1SJunchao Zhang for (i = 0; i < n; i++) {
6269566063dSJacob Faibussowitsch PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
6279371c9d4SSatish Balay if (!missing) {
6289371c9d4SSatish Balay *dups = PETSC_TRUE;
6299371c9d4SSatish Balay break;
6309371c9d4SSatish Balay }
631f1cab4e1SJunchao Zhang }
6329566063dSJacob Faibussowitsch PetscCall(PetscHSetIDestroy(&ht));
633f1cab4e1SJunchao Zhang }
6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
635f1cab4e1SJunchao Zhang }
636f1cab4e1SJunchao Zhang
637f1cab4e1SJunchao Zhang /*@
638811af0c4SBarry Smith PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`
639d2aeb606SJed Brown
640d2aeb606SJed Brown Not Collective
641d2aeb606SJed Brown
642d2aeb606SJed Brown Input Parameters:
643d2aeb606SJed Brown + key - the integer to locate
644d2aeb606SJed Brown . n - number of values in the array
6453334a463SJames Wright - X - array of `PetscMPIInt`
646d2aeb606SJed Brown
647d2aeb606SJed Brown Output Parameter:
648d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
649d2aeb606SJed Brown
650d2aeb606SJed Brown Level: intermediate
651d2aeb606SJed Brown
652db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
653d2aeb606SJed Brown @*/
PetscFindMPIInt(PetscMPIInt key,PetscCount n,const PetscMPIInt X[],PetscInt * loc)6546497c311SBarry Smith PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscCount n, const PetscMPIInt X[], PetscInt *loc)
655d71ae5a4SJacob Faibussowitsch {
6566497c311SBarry Smith PetscCount lo = 0, hi = n;
657d2aeb606SJed Brown
658d2aeb606SJed Brown PetscFunctionBegin;
6594f572ea9SToby Isaac PetscAssertPointer(loc, 4);
6609371c9d4SSatish Balay if (!n) {
6619371c9d4SSatish Balay *loc = -1;
6623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6639371c9d4SSatish Balay }
6644f572ea9SToby Isaac PetscAssertPointer(X, 3);
665d2aeb606SJed Brown while (hi - lo > 1) {
6666497c311SBarry Smith PetscCount mid = lo + (hi - lo) / 2;
66742f965a0SMatthew 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]);
66859e16d97SJunchao Zhang if (key < X[mid]) hi = mid;
669d2aeb606SJed Brown else lo = mid;
670d2aeb606SJed Brown }
671d76725a8SJose E. Roman PetscCall(PetscIntCast(key == X[lo] ? lo : -(lo + (MPIU_Count)(key > X[lo]) + 1), loc));
6723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
673e5c89e4eSSatish Balay }
674e5c89e4eSSatish Balay
675e5c89e4eSSatish Balay /*@
676811af0c4SBarry Smith PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
677811af0c4SBarry Smith changes a second array of `PetscInt` to match the sorted first array.
678e5c89e4eSSatish Balay
679e5c89e4eSSatish Balay Not Collective
680e5c89e4eSSatish Balay
681e5c89e4eSSatish Balay Input Parameters:
682e5c89e4eSSatish Balay + n - number of values
6833334a463SJames Wright . X - array of `PetscInt`
6843334a463SJames Wright - Y - second array of `PetscInt`
685e5c89e4eSSatish Balay
686e5c89e4eSSatish Balay Level: intermediate
687e5c89e4eSSatish Balay
688db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
689e5c89e4eSSatish Balay @*/
PetscSortIntWithArray(PetscCount n,PetscInt X[],PetscInt Y[])6906497c311SBarry Smith PetscErrorCode PetscSortIntWithArray(PetscCount n, PetscInt X[], PetscInt Y[])
691d71ae5a4SJacob Faibussowitsch {
69259e16d97SJunchao Zhang PetscInt pivot, t1, t2;
693e5c89e4eSSatish Balay
694e5c89e4eSSatish Balay PetscFunctionBegin;
6955f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
6963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
697c1f0200aSDmitry Karpeev }
698c1f0200aSDmitry Karpeev
699c1f0200aSDmitry Karpeev /*@
700811af0c4SBarry Smith PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
701811af0c4SBarry Smith changes a pair of `PetscInt` arrays to match the sorted first array.
702c1f0200aSDmitry Karpeev
703c1f0200aSDmitry Karpeev Not Collective
704c1f0200aSDmitry Karpeev
705c1f0200aSDmitry Karpeev Input Parameters:
706c1f0200aSDmitry Karpeev + n - number of values
7073334a463SJames Wright . X - array of `PestcInt`
7083334a463SJames Wright . Y - second array of `PestcInt` (first array of the pair)
7093334a463SJames Wright - Z - third array of `PestcInt` (second array of the pair)
710c1f0200aSDmitry Karpeev
711c1f0200aSDmitry Karpeev Level: intermediate
712c1f0200aSDmitry Karpeev
713db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
714c1f0200aSDmitry Karpeev @*/
PetscSortIntWithArrayPair(PetscCount n,PetscInt X[],PetscInt Y[],PetscInt Z[])7156497c311SBarry Smith PetscErrorCode PetscSortIntWithArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscInt Z[])
716d71ae5a4SJacob Faibussowitsch {
71759e16d97SJunchao Zhang PetscInt pivot, t1, t2, t3;
718c1f0200aSDmitry Karpeev
719c1f0200aSDmitry Karpeev PetscFunctionBegin;
7205f80ce2aSJacob Faibussowitsch QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
7213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
722c1f0200aSDmitry Karpeev }
723c1f0200aSDmitry Karpeev
72417d7d925SStefano Zampini /*@
7256497c311SBarry Smith PetscSortIntWithMPIIntArray - Sorts an array of `PetscInt` in place in increasing order;
7266497c311SBarry Smith changes a second array of `PetscMPI` to match the sorted first array.
7276497c311SBarry Smith
7286497c311SBarry Smith Not Collective
7296497c311SBarry Smith
7306497c311SBarry Smith Input Parameters:
7316497c311SBarry Smith + n - number of values
7323334a463SJames Wright . X - array of `PetscInt`
7333334a463SJames Wright - Y - second array of `PetscMPIInt`
7346497c311SBarry Smith
7356497c311SBarry Smith Level: intermediate
7366497c311SBarry Smith
73754c05997SPierre Jolivet .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
7386497c311SBarry Smith @*/
PetscSortIntWithMPIIntArray(PetscCount n,PetscInt X[],PetscMPIInt Y[])7396497c311SBarry Smith PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount n, PetscInt X[], PetscMPIInt Y[])
7406497c311SBarry Smith {
7416497c311SBarry Smith PetscInt pivot, t1;
7426497c311SBarry Smith PetscMPIInt t2;
7436497c311SBarry Smith
7446497c311SBarry Smith PetscFunctionBegin;
7456497c311SBarry Smith QuickSort2(PetscSortIntWithMPIIntArray, X, Y, n, pivot, t1, t2);
7466497c311SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
7476497c311SBarry Smith }
7486497c311SBarry Smith
7496497c311SBarry Smith /*@
750811af0c4SBarry Smith PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
751811af0c4SBarry Smith changes a second array of `PetscCount` to match the sorted first array.
752981bb840SJunchao Zhang
753981bb840SJunchao Zhang Not Collective
754981bb840SJunchao Zhang
755981bb840SJunchao Zhang Input Parameters:
756981bb840SJunchao Zhang + n - number of values
7573334a463SJames Wright . X - array of `PetscInt`
7583334a463SJames Wright - Y - second array of `PetscCount`
759981bb840SJunchao Zhang
760981bb840SJunchao Zhang Level: intermediate
761981bb840SJunchao Zhang
76254c05997SPierre Jolivet .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
763981bb840SJunchao Zhang @*/
PetscSortIntWithCountArray(PetscCount n,PetscInt X[],PetscCount Y[])764d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
765d71ae5a4SJacob Faibussowitsch {
766981bb840SJunchao Zhang PetscInt pivot, t1;
767981bb840SJunchao Zhang PetscCount t2;
768981bb840SJunchao Zhang
769981bb840SJunchao Zhang PetscFunctionBegin;
7705f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
7713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
772981bb840SJunchao Zhang }
773981bb840SJunchao Zhang
774981bb840SJunchao Zhang /*@
775811af0c4SBarry Smith PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
776811af0c4SBarry Smith changes a `PetscInt` array and a `PetscCount` array to match the sorted first array.
777981bb840SJunchao Zhang
778981bb840SJunchao Zhang Not Collective
779981bb840SJunchao Zhang
780981bb840SJunchao Zhang Input Parameters:
781981bb840SJunchao Zhang + n - number of values
7823334a463SJames Wright . X - array of `PetscInt`
7833334a463SJames Wright . Y - second array of `PetscInt` (first array of the pair)
7843334a463SJames Wright - Z - third array of `PetscCount` (second array of the pair)
785981bb840SJunchao Zhang
786981bb840SJunchao Zhang Level: intermediate
787981bb840SJunchao Zhang
788811af0c4SBarry Smith Note:
789981bb840SJunchao 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.
790981bb840SJunchao Zhang
79154c05997SPierre Jolivet .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
792981bb840SJunchao Zhang @*/
PetscSortIntWithIntCountArrayPair(PetscCount n,PetscInt X[],PetscInt Y[],PetscCount Z[])793d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
794d71ae5a4SJacob Faibussowitsch {
795981bb840SJunchao Zhang PetscInt pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
796981bb840SJunchao Zhang PetscCount t3; /* temp for Z[] */
797981bb840SJunchao Zhang
798981bb840SJunchao Zhang PetscFunctionBegin;
7995f80ce2aSJacob Faibussowitsch QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
8003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
801981bb840SJunchao Zhang }
802981bb840SJunchao Zhang
803981bb840SJunchao Zhang /*@
804811af0c4SBarry Smith PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.
8056a7c706bSVaclav Hapla
8066a7c706bSVaclav Hapla Not Collective
8076a7c706bSVaclav Hapla
8086a7c706bSVaclav Hapla Input Parameters:
8096a7c706bSVaclav Hapla + n - number of values
8103334a463SJames Wright - X - array of `PetscMPIInt`
8116a7c706bSVaclav Hapla
8122fe279fdSBarry Smith Output Parameter:
8136a7c706bSVaclav Hapla . sorted - flag whether the array is sorted
8146a7c706bSVaclav Hapla
8156a7c706bSVaclav Hapla Level: intermediate
8166a7c706bSVaclav Hapla
817db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
8186a7c706bSVaclav Hapla @*/
PetscSortedMPIInt(PetscCount n,const PetscMPIInt X[],PetscBool * sorted)8196497c311SBarry Smith PetscErrorCode PetscSortedMPIInt(PetscCount n, const PetscMPIInt X[], PetscBool *sorted)
820d71ae5a4SJacob Faibussowitsch {
8216a7c706bSVaclav Hapla PetscFunctionBegin;
8226a7c706bSVaclav Hapla PetscSorted(n, X, *sorted);
8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8246a7c706bSVaclav Hapla }
8256a7c706bSVaclav Hapla
8266a7c706bSVaclav Hapla /*@
827811af0c4SBarry Smith PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.
82817d7d925SStefano Zampini
82917d7d925SStefano Zampini Not Collective
83017d7d925SStefano Zampini
83117d7d925SStefano Zampini Input Parameters:
83217d7d925SStefano Zampini + n - number of values
8333334a463SJames Wright - X - array of `PetscMPIInt`
83417d7d925SStefano Zampini
83517d7d925SStefano Zampini Level: intermediate
83617d7d925SStefano Zampini
837811af0c4SBarry Smith Note:
838676f2a66SJacob Faibussowitsch This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
839a5b23f4aSJose E. Roman is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
840676f2a66SJacob Faibussowitsch code to see which routine is fastest.
841676f2a66SJacob Faibussowitsch
842db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
84317d7d925SStefano Zampini @*/
PetscSortMPIInt(PetscCount n,PetscMPIInt X[])8446497c311SBarry Smith PetscErrorCode PetscSortMPIInt(PetscCount n, PetscMPIInt X[])
845d71ae5a4SJacob Faibussowitsch {
84659e16d97SJunchao Zhang PetscMPIInt pivot, t1;
84717d7d925SStefano Zampini
84817d7d925SStefano Zampini PetscFunctionBegin;
8495f80ce2aSJacob Faibussowitsch QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
8503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
85117d7d925SStefano Zampini }
85217d7d925SStefano Zampini
85317d7d925SStefano Zampini /*@
854811af0c4SBarry Smith PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries
85517d7d925SStefano Zampini
85617d7d925SStefano Zampini Not Collective
85717d7d925SStefano Zampini
85817d7d925SStefano Zampini Input Parameters:
85917d7d925SStefano Zampini + n - number of values
8603334a463SJames Wright - X - array of `PetscMPIInt`
86117d7d925SStefano Zampini
86217d7d925SStefano Zampini Output Parameter:
86317d7d925SStefano Zampini . n - number of non-redundant values
86417d7d925SStefano Zampini
86517d7d925SStefano Zampini Level: intermediate
86617d7d925SStefano Zampini
867db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
86817d7d925SStefano Zampini @*/
PetscSortRemoveDupsMPIInt(PetscInt * n,PetscMPIInt X[])869d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
870d71ae5a4SJacob Faibussowitsch {
8715f80ce2aSJacob Faibussowitsch PetscInt s = 0, N = *n, b = 0;
87217d7d925SStefano Zampini
87317d7d925SStefano Zampini PetscFunctionBegin;
8749566063dSJacob Faibussowitsch PetscCall(PetscSortMPIInt(N, X));
8755f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < N - 1; i++) {
87659e16d97SJunchao Zhang if (X[b + s + 1] != X[b]) {
8779371c9d4SSatish Balay X[b + 1] = X[b + s + 1];
8789371c9d4SSatish Balay b++;
879a297a907SKarl Rupp } else s++;
88017d7d925SStefano Zampini }
88117d7d925SStefano Zampini *n = N - s;
8823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
88317d7d925SStefano Zampini }
884c1f0200aSDmitry Karpeev
8854d615ea0SBarry Smith /*@
886811af0c4SBarry Smith PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
887811af0c4SBarry Smith changes a second `PetscMPIInt` array to match the sorted first array.
8884d615ea0SBarry Smith
8894d615ea0SBarry Smith Not Collective
8904d615ea0SBarry Smith
8914d615ea0SBarry Smith Input Parameters:
8924d615ea0SBarry Smith + n - number of values
8933334a463SJames Wright . X - array of `PetscMPIInt`
8943334a463SJames Wright - Y - second array of `PetscMPIInt`
8954d615ea0SBarry Smith
8964d615ea0SBarry Smith Level: intermediate
8974d615ea0SBarry Smith
898db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
8994d615ea0SBarry Smith @*/
PetscSortMPIIntWithArray(PetscCount n,PetscMPIInt X[],PetscMPIInt Y[])9006497c311SBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscCount n, PetscMPIInt X[], PetscMPIInt Y[])
901d71ae5a4SJacob Faibussowitsch {
90259e16d97SJunchao Zhang PetscMPIInt pivot, t1, t2;
9034d615ea0SBarry Smith
9044d615ea0SBarry Smith PetscFunctionBegin;
9055f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
9063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9075569a785SJunchao Zhang }
9085569a785SJunchao Zhang
9095569a785SJunchao Zhang /*@
910811af0c4SBarry Smith PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
911811af0c4SBarry Smith changes a second array of `PetscInt` to match the sorted first array.
9125569a785SJunchao Zhang
9135569a785SJunchao Zhang Not Collective
9145569a785SJunchao Zhang
9155569a785SJunchao Zhang Input Parameters:
9165569a785SJunchao Zhang + n - number of values
9173334a463SJames Wright . X - array of `PetscMPIInt`
9183334a463SJames Wright - Y - second array of `PetscInt`
9195569a785SJunchao Zhang
9205569a785SJunchao Zhang Level: intermediate
9215569a785SJunchao Zhang
922811af0c4SBarry Smith Note:
923811af0c4SBarry Smith This routine is useful when one needs to sort MPI ranks with other integer arrays.
9245569a785SJunchao Zhang
925db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
9265569a785SJunchao Zhang @*/
PetscSortMPIIntWithIntArray(PetscCount n,PetscMPIInt X[],PetscInt Y[])9276497c311SBarry Smith PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount n, PetscMPIInt X[], PetscInt Y[])
928d71ae5a4SJacob Faibussowitsch {
92959e16d97SJunchao Zhang PetscMPIInt pivot, t1;
9305569a785SJunchao Zhang PetscInt t2;
9315569a785SJunchao Zhang
9325569a785SJunchao Zhang PetscFunctionBegin;
9335f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
9343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
935e5c89e4eSSatish Balay }
936e5c89e4eSSatish Balay
937e5c89e4eSSatish Balay /*@
938811af0c4SBarry Smith PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
939811af0c4SBarry Smith changes a second `PetscScalar` array to match the sorted first array.
940e5c89e4eSSatish Balay
941e5c89e4eSSatish Balay Not Collective
942e5c89e4eSSatish Balay
943e5c89e4eSSatish Balay Input Parameters:
944e5c89e4eSSatish Balay + n - number of values
9453334a463SJames Wright . X - array of `PetscInt`
9463334a463SJames Wright - Y - second array of `PetscScalar`
947e5c89e4eSSatish Balay
948e5c89e4eSSatish Balay Level: intermediate
949e5c89e4eSSatish Balay
950db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
951e5c89e4eSSatish Balay @*/
PetscSortIntWithScalarArray(PetscCount n,PetscInt X[],PetscScalar Y[])9526497c311SBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscCount n, PetscInt X[], PetscScalar Y[])
953d71ae5a4SJacob Faibussowitsch {
95459e16d97SJunchao Zhang PetscInt pivot, t1;
95559e16d97SJunchao Zhang PetscScalar t2;
956e5c89e4eSSatish Balay
957e5c89e4eSSatish Balay PetscFunctionBegin;
9585f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
9593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
96017df18f2SToby Isaac }
96117df18f2SToby Isaac
9626c2863d0SBarry Smith /*@C
963811af0c4SBarry Smith PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
96417df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must
96517df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting.
96617df18f2SToby Isaac
967cc4c1da9SBarry Smith Not Collective, No Fortran Support
96817df18f2SToby Isaac
96917df18f2SToby Isaac Input Parameters:
97017df18f2SToby Isaac + n - number of values
9713334a463SJames Wright . X - array of `PetscInt`
97259e16d97SJunchao Zhang . Y - second array of data
97317df18f2SToby Isaac . size - sizeof elements in the data array in bytes
97459e16d97SJunchao Zhang - t2 - workspace of "size" bytes used when sorting
97517df18f2SToby Isaac
97617df18f2SToby Isaac Level: intermediate
97717df18f2SToby Isaac
978db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
97917df18f2SToby Isaac @*/
PetscSortIntWithDataArray(PetscCount n,PetscInt X[],void * Y,size_t size,void * t2)9806497c311SBarry Smith PetscErrorCode PetscSortIntWithDataArray(PetscCount n, PetscInt X[], void *Y, size_t size, void *t2)
981d71ae5a4SJacob Faibussowitsch {
98259e16d97SJunchao Zhang char *YY = (char *)Y;
9836497c311SBarry Smith PetscCount hi = n - 1;
9846497c311SBarry Smith PetscInt pivot, t1;
98517df18f2SToby Isaac
98617df18f2SToby Isaac PetscFunctionBegin;
98717df18f2SToby Isaac if (n < 8) {
9886497c311SBarry Smith for (PetscCount i = 0; i < n; i++) {
98959e16d97SJunchao Zhang pivot = X[i];
9906497c311SBarry Smith for (PetscCount j = i + 1; j < n; j++) {
99159e16d97SJunchao Zhang if (pivot > X[j]) {
99259e16d97SJunchao Zhang SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
99359e16d97SJunchao Zhang pivot = X[i];
99417df18f2SToby Isaac }
99517df18f2SToby Isaac }
99617df18f2SToby Isaac }
99717df18f2SToby Isaac } else {
99859e16d97SJunchao Zhang /* Two way partition */
9996497c311SBarry Smith PetscCount l = 0, r = hi;
10005f80ce2aSJacob Faibussowitsch
10015f80ce2aSJacob Faibussowitsch pivot = X[MEDIAN(X, hi)];
100259e16d97SJunchao Zhang while (1) {
100359e16d97SJunchao Zhang while (X[l] < pivot) l++;
100459e16d97SJunchao Zhang while (X[r] > pivot) r--;
10059371c9d4SSatish Balay if (l >= r) {
10069371c9d4SSatish Balay r++;
10079371c9d4SSatish Balay break;
10089371c9d4SSatish Balay }
100959e16d97SJunchao Zhang SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
101059e16d97SJunchao Zhang l++;
101159e16d97SJunchao Zhang r--;
101259e16d97SJunchao Zhang }
10139566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
10149566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
101517df18f2SToby Isaac }
10163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
101717df18f2SToby Isaac }
101817df18f2SToby Isaac
101921e72a00SBarry Smith /*@
1020811af0c4SBarry Smith PetscMergeIntArray - Merges two SORTED `PetscInt` arrays, removes duplicate elements.
102121e72a00SBarry Smith
102221e72a00SBarry Smith Not Collective
102321e72a00SBarry Smith
102421e72a00SBarry Smith Input Parameters:
102521e72a00SBarry Smith + an - number of values in the first array
10263334a463SJames Wright . aI - first sorted array of `PetscInt`
102721e72a00SBarry Smith . bn - number of values in the second array
10283334a463SJames Wright - bI - second array of `PetscInt`
102921e72a00SBarry Smith
103021e72a00SBarry Smith Output Parameters:
103121e72a00SBarry Smith + n - number of values in the merged array
10326c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided
103321e72a00SBarry Smith
103421e72a00SBarry Smith Level: intermediate
103521e72a00SBarry Smith
1036db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
103721e72a00SBarry Smith @*/
PetscMergeIntArray(PetscInt an,const PetscInt aI[],PetscInt bn,const PetscInt bI[],PetscInt * n,PetscInt ** L)1038d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
1039d71ae5a4SJacob Faibussowitsch {
104021e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k;
104121e72a00SBarry Smith
1042362febeeSStefano Zampini PetscFunctionBegin;
104321e72a00SBarry Smith if (!L_) {
10449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(an + bn, L));
104521e72a00SBarry Smith L_ = *L;
104621e72a00SBarry Smith }
104721e72a00SBarry Smith k = ak = bk = 0;
104821e72a00SBarry Smith while (ak < an && bk < bn) {
104921e72a00SBarry Smith if (aI[ak] == bI[bk]) {
105021e72a00SBarry Smith L_[k] = aI[ak];
105121e72a00SBarry Smith ++ak;
105221e72a00SBarry Smith ++bk;
105321e72a00SBarry Smith ++k;
105421e72a00SBarry Smith } else if (aI[ak] < bI[bk]) {
105521e72a00SBarry Smith L_[k] = aI[ak];
105621e72a00SBarry Smith ++ak;
105721e72a00SBarry Smith ++k;
105821e72a00SBarry Smith } else {
105921e72a00SBarry Smith L_[k] = bI[bk];
106021e72a00SBarry Smith ++bk;
106121e72a00SBarry Smith ++k;
106221e72a00SBarry Smith }
106321e72a00SBarry Smith }
106421e72a00SBarry Smith if (ak < an) {
10659566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
106621e72a00SBarry Smith k += (an - ak);
106721e72a00SBarry Smith }
106821e72a00SBarry Smith if (bk < bn) {
10699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
107021e72a00SBarry Smith k += (bn - bk);
107121e72a00SBarry Smith }
107221e72a00SBarry Smith *n = k;
10733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
107421e72a00SBarry Smith }
1075b4301105SBarry Smith
1076c1f0200aSDmitry Karpeev /*@
1077811af0c4SBarry Smith PetscMergeIntArrayPair - Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
1078c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged
1079c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair.
1080c1f0200aSDmitry Karpeev
1081c1f0200aSDmitry Karpeev Not Collective
1082c1f0200aSDmitry Karpeev
1083c1f0200aSDmitry Karpeev Input Parameters:
1084c1f0200aSDmitry Karpeev + an - number of values in the first array
10853334a463SJames Wright . aI - first sorted array of `PetscInt`
10863334a463SJames Wright . aJ - first additional array of `PetscInt`
1087c1f0200aSDmitry Karpeev . bn - number of values in the second array
10883334a463SJames Wright . bI - second array of `PetscInt`
10893334a463SJames Wright - bJ - second additional of `PetscInt`
1090c1f0200aSDmitry Karpeev
1091c1f0200aSDmitry Karpeev Output Parameters:
1092c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn)
109314c49006SJed Brown . L - merged sorted array
1094c1f0200aSDmitry Karpeev - J - merged additional array
1095c1f0200aSDmitry Karpeev
1096811af0c4SBarry Smith Note:
1097da81f932SPierre 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
1098811af0c4SBarry Smith
1099c1f0200aSDmitry Karpeev Level: intermediate
1100c1f0200aSDmitry Karpeev
1101db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1102c1f0200aSDmitry Karpeev @*/
PetscMergeIntArrayPair(PetscInt an,const PetscInt aI[],const PetscInt aJ[],PetscInt bn,const PetscInt bI[],const PetscInt bJ[],PetscInt * n,PetscInt * L[],PetscInt * J[])1103ce78bad3SBarry Smith PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt *L[], PetscInt *J[])
1104d71ae5a4SJacob Faibussowitsch {
110570d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k;
1106c1f0200aSDmitry Karpeev
110714c49006SJed Brown PetscFunctionBegin;
11084f572ea9SToby Isaac PetscAssertPointer(L, 8);
11094f572ea9SToby Isaac PetscAssertPointer(J, 9);
1110c1f0200aSDmitry Karpeev n_ = an + bn;
1111c1f0200aSDmitry Karpeev *n = n_;
111248a46eb9SPierre Jolivet if (!*L) PetscCall(PetscMalloc1(n_, L));
1113d7aa01a8SHong Zhang L_ = *L;
111448a46eb9SPierre Jolivet if (!*J) PetscCall(PetscMalloc1(n_, J));
1115c1f0200aSDmitry Karpeev J_ = *J;
1116c1f0200aSDmitry Karpeev k = ak = bk = 0;
1117c1f0200aSDmitry Karpeev while (ak < an && bk < bn) {
1118c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) {
1119d7aa01a8SHong Zhang L_[k] = aI[ak];
1120c1f0200aSDmitry Karpeev J_[k] = aJ[ak];
1121c1f0200aSDmitry Karpeev ++ak;
1122c1f0200aSDmitry Karpeev ++k;
1123a297a907SKarl Rupp } else {
1124d7aa01a8SHong Zhang L_[k] = bI[bk];
1125c1f0200aSDmitry Karpeev J_[k] = bJ[bk];
1126c1f0200aSDmitry Karpeev ++bk;
1127c1f0200aSDmitry Karpeev ++k;
1128c1f0200aSDmitry Karpeev }
1129c1f0200aSDmitry Karpeev }
1130c1f0200aSDmitry Karpeev if (ak < an) {
11319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
11329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1133c1f0200aSDmitry Karpeev k += (an - ak);
1134c1f0200aSDmitry Karpeev }
1135c1f0200aSDmitry Karpeev if (bk < bn) {
11369566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
11379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1138c1f0200aSDmitry Karpeev }
11393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1140c1f0200aSDmitry Karpeev }
1141c1f0200aSDmitry Karpeev
1142e498c390SJed Brown /*@
1143811af0c4SBarry Smith PetscMergeMPIIntArray - Merges two SORTED `PetscMPIInt` arrays.
1144e498c390SJed Brown
1145e498c390SJed Brown Not Collective
1146e498c390SJed Brown
1147e498c390SJed Brown Input Parameters:
1148e498c390SJed Brown + an - number of values in the first array
11493334a463SJames Wright . aI - first sorted array of `PetscMPIInt`
1150e498c390SJed Brown . bn - number of values in the second array
11513334a463SJames Wright - bI - second array of `PetscMPIInt`
1152e498c390SJed Brown
1153e498c390SJed Brown Output Parameters:
1154e498c390SJed Brown + n - number of values in the merged array (<= an + bn)
1155e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed
1156e498c390SJed Brown
1157e498c390SJed Brown Level: intermediate
1158e498c390SJed Brown
1159db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1160e498c390SJed Brown @*/
PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt * n,PetscMPIInt ** L)1161d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1162d71ae5a4SJacob Faibussowitsch {
1163e498c390SJed Brown PetscInt ai, bi, k;
1164e498c390SJed Brown
1165e498c390SJed Brown PetscFunctionBegin;
116657508eceSPierre Jolivet if (!*L) PetscCall(PetscMalloc1(an + bn, L));
1167e498c390SJed Brown for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1168e498c390SJed Brown PetscInt t = -1;
1169e498c390SJed Brown for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
1170fbccb6d4SPierre Jolivet for (; bi < bn && bI[bi] == t; bi++);
1171e498c390SJed Brown for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
1172fbccb6d4SPierre Jolivet for (; ai < an && aI[ai] == t; ai++);
1173e498c390SJed Brown }
1174e498c390SJed Brown *n = k;
11753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1176e498c390SJed Brown }
1177e5c89e4eSSatish Balay
11786c2863d0SBarry Smith /*@C
117942eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically
118042eaa7f3SBarry Smith
1181cc4c1da9SBarry Smith Not Collective, No Fortran Support
118242eaa7f3SBarry Smith
118342eaa7f3SBarry Smith Input Parameters:
118442eaa7f3SBarry Smith + n - number of values
118542eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked
118642eaa7f3SBarry Smith - parentid - indicates the parent of each entry
118742eaa7f3SBarry Smith
118842eaa7f3SBarry Smith Output Parameters:
118942eaa7f3SBarry Smith + Nlevels - the number of levels
119042eaa7f3SBarry Smith . Level - for each node tells its level
1191aec76313SJacob Faibussowitsch . Levelcnt - the number of nodes on each level
119242eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc
119342eaa7f3SBarry Smith - Column - for each id tells its column index
119442eaa7f3SBarry Smith
11956c2863d0SBarry Smith Level: developer
119642eaa7f3SBarry Smith
1197811af0c4SBarry Smith Note:
119895452b02SPatrick Sanan This code is not currently used
119942eaa7f3SBarry Smith
1200db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
120142eaa7f3SBarry Smith @*/
PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt * Nlevels,PetscInt ** Level,PetscInt ** Levelcnt,PetscInt ** Idbylevel,PetscInt ** Column)1202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1203d71ae5a4SJacob Faibussowitsch {
120442eaa7f3SBarry Smith PetscInt i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1205ace3abfcSBarry Smith PetscBool done = PETSC_FALSE;
120642eaa7f3SBarry Smith
120742eaa7f3SBarry Smith PetscFunctionBegin;
120808401ef6SPierre Jolivet PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
120942eaa7f3SBarry Smith for (i = 0; i < n; i++) {
121042eaa7f3SBarry Smith if (mask[i]) continue;
121108401ef6SPierre Jolivet PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1212cc73adaaSBarry Smith PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
121342eaa7f3SBarry Smith }
121442eaa7f3SBarry Smith
121542eaa7f3SBarry Smith for (i = 0; i < n; i++) {
121642eaa7f3SBarry Smith if (!mask[i]) nmask++;
121742eaa7f3SBarry Smith }
121842eaa7f3SBarry Smith
121942eaa7f3SBarry Smith /* determine the level in the tree of each node */
12209566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n, &level));
1221a297a907SKarl Rupp
122242eaa7f3SBarry Smith level[0] = 1;
122342eaa7f3SBarry Smith while (!done) {
122442eaa7f3SBarry Smith done = PETSC_TRUE;
122542eaa7f3SBarry Smith for (i = 0; i < n; i++) {
122642eaa7f3SBarry Smith if (mask[i]) continue;
122742eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
122842eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE;
122942eaa7f3SBarry Smith }
123042eaa7f3SBarry Smith }
123142eaa7f3SBarry Smith for (i = 0; i < n; i++) {
123242eaa7f3SBarry Smith level[i]--;
123342eaa7f3SBarry Smith nlevels = PetscMax(nlevels, level[i]);
123442eaa7f3SBarry Smith }
123542eaa7f3SBarry Smith
123642eaa7f3SBarry Smith /* count the number of nodes on each level and its max */
12379566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nlevels, &levelcnt));
123842eaa7f3SBarry Smith for (i = 0; i < n; i++) {
123942eaa7f3SBarry Smith if (mask[i]) continue;
124042eaa7f3SBarry Smith levelcnt[level[i] - 1]++;
124142eaa7f3SBarry Smith }
1242a297a907SKarl Rupp for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);
124342eaa7f3SBarry Smith
124442eaa7f3SBarry Smith /* for each level sort the ids by the parent id */
12459566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
12469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nmask, &idbylevel));
124742eaa7f3SBarry Smith for (j = 1; j <= nlevels; j++) {
124842eaa7f3SBarry Smith cnt = 0;
124942eaa7f3SBarry Smith for (i = 0; i < n; i++) {
125042eaa7f3SBarry Smith if (mask[i]) continue;
125142eaa7f3SBarry Smith if (level[i] != j) continue;
125242eaa7f3SBarry Smith workid[cnt] = i;
125342eaa7f3SBarry Smith workparentid[cnt++] = parentid[i];
125442eaa7f3SBarry Smith }
125542eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0);
125642eaa7f3SBarry Smith PetscIntView(cnt,workid,0);
12579566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
125842eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0);
125942eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/
12609566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
126142eaa7f3SBarry Smith tcnt += cnt;
126242eaa7f3SBarry Smith }
126308401ef6SPierre Jolivet PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
12649566063dSJacob Faibussowitsch PetscCall(PetscFree2(workid, workparentid));
126542eaa7f3SBarry Smith
126642eaa7f3SBarry Smith /* for each node list its column */
12679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &column));
126842eaa7f3SBarry Smith cnt = 0;
126942eaa7f3SBarry Smith for (j = 0; j < nlevels; j++) {
1270ad540459SPierre Jolivet for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
127142eaa7f3SBarry Smith }
127242eaa7f3SBarry Smith
127342eaa7f3SBarry Smith *Nlevels = nlevels;
127442eaa7f3SBarry Smith *Level = level;
127542eaa7f3SBarry Smith *Levelcnt = levelcnt;
127642eaa7f3SBarry Smith *Idbylevel = idbylevel;
127742eaa7f3SBarry Smith *Column = column;
12783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
127942eaa7f3SBarry Smith }
1280ce605777SToby Isaac
1281ce605777SToby Isaac /*@
1282811af0c4SBarry Smith PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.
1283ce605777SToby Isaac
1284ce605777SToby Isaac Collective
1285ce605777SToby Isaac
1286ce605777SToby Isaac Input Parameters:
1287ce605777SToby Isaac + comm - the MPI communicator
12883334a463SJames Wright . n - the local number of `PetscInt`
12893334a463SJames Wright - keys - the local array of `PetscInt`
1290ce605777SToby Isaac
1291ce605777SToby Isaac Output Parameters:
1292ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1293ce605777SToby Isaac
1294ce605777SToby Isaac Level: developer
1295ce605777SToby Isaac
1296db781477SPatrick Sanan .seealso: `PetscParallelSortInt()`
1297ce605777SToby Isaac @*/
PetscParallelSortedInt(MPI_Comm comm,PetscInt n,const PetscInt keys[],PetscBool * is_sorted)1298d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1299d71ae5a4SJacob Faibussowitsch {
1300ce605777SToby Isaac PetscBool sorted;
1301ce605777SToby Isaac PetscInt i, min, max, prevmax;
13022a1da528SToby Isaac PetscMPIInt rank;
1303ce605777SToby Isaac
1304ce605777SToby Isaac PetscFunctionBegin;
1305ce605777SToby Isaac sorted = PETSC_TRUE;
13061690c2aeSBarry Smith min = PETSC_INT_MAX;
13071690c2aeSBarry Smith max = PETSC_INT_MIN;
1308ce605777SToby Isaac if (n) {
1309ce605777SToby Isaac min = keys[0];
1310ce605777SToby Isaac max = keys[0];
1311ce605777SToby Isaac }
1312ce605777SToby Isaac for (i = 1; i < n; i++) {
1313ce605777SToby Isaac if (keys[i] < keys[i - 1]) break;
1314ce605777SToby Isaac min = PetscMin(min, keys[i]);
1315ce605777SToby Isaac max = PetscMax(max, keys[i]);
1316ce605777SToby Isaac }
1317ce605777SToby Isaac if (i < n) sorted = PETSC_FALSE;
13181690c2aeSBarry Smith prevmax = PETSC_INT_MIN;
13199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
13209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
13211690c2aeSBarry Smith if (rank == 0) prevmax = PETSC_INT_MIN;
1322ce605777SToby Isaac if (prevmax > min) sorted = PETSC_FALSE;
1323*5440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&sorted, is_sorted, 1, MPI_C_BOOL, MPI_LAND, comm));
13243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1325ce605777SToby Isaac }
1326