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 65*811af0c4SBarry 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 98*811af0c4SBarry 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 /*@ 244*811af0c4SBarry 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 @*/ 2599371c9d4SSatish Balay PetscErrorCode PetscSortedInt(PetscInt n, const PetscInt X[], PetscBool *sorted) { 2606a7c706bSVaclav Hapla PetscFunctionBegin; 2615f80ce2aSJacob Faibussowitsch if (n) PetscValidIntPointer(X, 2); 2625f80ce2aSJacob Faibussowitsch PetscValidBoolPointer(sorted, 3); 2636a7c706bSVaclav Hapla PetscSorted(n, X, *sorted); 2646a7c706bSVaclav Hapla PetscFunctionReturn(0); 2656a7c706bSVaclav Hapla } 2666a7c706bSVaclav Hapla 2676a7c706bSVaclav Hapla /*@ 268*811af0c4SBarry Smith PetscSortInt - Sorts an array of `PetscInt` in place in increasing order. 269e5c89e4eSSatish Balay 270e5c89e4eSSatish Balay Not Collective 271e5c89e4eSSatish Balay 272e5c89e4eSSatish Balay Input Parameters: 273e5c89e4eSSatish Balay + n - number of values 27459e16d97SJunchao Zhang - X - array of integers 275e5c89e4eSSatish Balay 276*811af0c4SBarry Smith Note: 277*811af0c4SBarry Smith This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array 278a5b23f4aSJose E. Roman is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their 279676f2a66SJacob Faibussowitsch code to see which routine is fastest. 280676f2a66SJacob Faibussowitsch 281e5c89e4eSSatish Balay Level: intermediate 282e5c89e4eSSatish Balay 283db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()` 284e5c89e4eSSatish Balay @*/ 2859371c9d4SSatish Balay PetscErrorCode PetscSortInt(PetscInt n, PetscInt X[]) { 28659e16d97SJunchao Zhang PetscInt pivot, t1; 287e5c89e4eSSatish Balay 288e5c89e4eSSatish Balay PetscFunctionBegin; 2895f80ce2aSJacob Faibussowitsch if (n) PetscValidIntPointer(X, 2); 2905f80ce2aSJacob Faibussowitsch QuickSort1(PetscSortInt, X, n, pivot, t1); 291e5c89e4eSSatish Balay PetscFunctionReturn(0); 292e5c89e4eSSatish Balay } 293e5c89e4eSSatish Balay 2941db36a52SBarry Smith /*@ 295*811af0c4SBarry Smith PetscSortReverseInt - Sorts an array of `PetscInt` in place in decreasing order. 296ce605777SToby Isaac 297ce605777SToby Isaac Not Collective 298ce605777SToby Isaac 299ce605777SToby Isaac Input Parameters: 300ce605777SToby Isaac + n - number of values 301ce605777SToby Isaac - X - array of integers 302ce605777SToby Isaac 303ce605777SToby Isaac Level: intermediate 304ce605777SToby Isaac 305db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()` 306ce605777SToby Isaac @*/ 3079371c9d4SSatish Balay PetscErrorCode PetscSortReverseInt(PetscInt n, PetscInt X[]) { 308ce605777SToby Isaac PetscInt pivot, t1; 309ce605777SToby Isaac 310ce605777SToby Isaac PetscFunctionBegin; 3115f80ce2aSJacob Faibussowitsch if (n) PetscValidIntPointer(X, 2); 3125f80ce2aSJacob Faibussowitsch QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1); 313ce605777SToby Isaac PetscFunctionReturn(0); 314ce605777SToby Isaac } 315ce605777SToby Isaac 316ce605777SToby Isaac /*@ 317*811af0c4SBarry Smith PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array 31822ab5688SLisandro Dalcin 31922ab5688SLisandro Dalcin Not Collective 32022ab5688SLisandro Dalcin 32122ab5688SLisandro Dalcin Input Parameters: 32222ab5688SLisandro Dalcin + n - number of values 32359e16d97SJunchao Zhang - X - sorted array of integers 32422ab5688SLisandro Dalcin 32522ab5688SLisandro Dalcin Output Parameter: 32622ab5688SLisandro Dalcin . n - number of non-redundant values 32722ab5688SLisandro Dalcin 32822ab5688SLisandro Dalcin Level: intermediate 32922ab5688SLisandro Dalcin 330db781477SPatrick Sanan .seealso: `PetscSortInt()` 33122ab5688SLisandro Dalcin @*/ 3329371c9d4SSatish Balay PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[]) { 33322ab5688SLisandro Dalcin PetscInt i, s = 0, N = *n, b = 0; 33422ab5688SLisandro Dalcin 33522ab5688SLisandro Dalcin PetscFunctionBegin; 3365f80ce2aSJacob Faibussowitsch PetscValidIntPointer(n, 1); 337fb61b9e4SStefano Zampini PetscCheckSorted(*n, X); 33822ab5688SLisandro Dalcin for (i = 0; i < N - 1; i++) { 33959e16d97SJunchao Zhang if (X[b + s + 1] != X[b]) { 3409371c9d4SSatish Balay X[b + 1] = X[b + s + 1]; 3419371c9d4SSatish Balay b++; 34222ab5688SLisandro Dalcin } else s++; 34322ab5688SLisandro Dalcin } 34422ab5688SLisandro Dalcin *n = N - s; 34522ab5688SLisandro Dalcin PetscFunctionReturn(0); 34622ab5688SLisandro Dalcin } 34722ab5688SLisandro Dalcin 34822ab5688SLisandro Dalcin /*@ 349*811af0c4SBarry Smith PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates 350b6a18d21SVaclav Hapla 351b6a18d21SVaclav Hapla Not Collective 352b6a18d21SVaclav Hapla 353b6a18d21SVaclav Hapla Input Parameters: 354b6a18d21SVaclav Hapla + n - number of values 355b6a18d21SVaclav Hapla - X - sorted array of integers 356b6a18d21SVaclav Hapla 357b6a18d21SVaclav Hapla Output Parameter: 358b6a18d21SVaclav Hapla . dups - True if the array has dups, otherwise false 359b6a18d21SVaclav Hapla 360b6a18d21SVaclav Hapla Level: intermediate 361b6a18d21SVaclav Hapla 362db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()` 363b6a18d21SVaclav Hapla @*/ 3649371c9d4SSatish Balay PetscErrorCode PetscSortedCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *flg) { 365b6a18d21SVaclav Hapla PetscInt i; 366b6a18d21SVaclav Hapla 367b6a18d21SVaclav Hapla PetscFunctionBegin; 368b6a18d21SVaclav Hapla PetscCheckSorted(n, X); 369b6a18d21SVaclav Hapla *flg = PETSC_FALSE; 370b6a18d21SVaclav Hapla for (i = 0; i < n - 1; i++) { 371b6a18d21SVaclav Hapla if (X[i + 1] == X[i]) { 372b6a18d21SVaclav Hapla *flg = PETSC_TRUE; 373b6a18d21SVaclav Hapla break; 374b6a18d21SVaclav Hapla } 375b6a18d21SVaclav Hapla } 376b6a18d21SVaclav Hapla PetscFunctionReturn(0); 377b6a18d21SVaclav Hapla } 378b6a18d21SVaclav Hapla 379b6a18d21SVaclav Hapla /*@ 380*811af0c4SBarry Smith PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries 3811db36a52SBarry Smith 3821db36a52SBarry Smith Not Collective 3831db36a52SBarry Smith 3841db36a52SBarry Smith Input Parameters: 3851db36a52SBarry Smith + n - number of values 38659e16d97SJunchao Zhang - X - array of integers 3871db36a52SBarry Smith 3881db36a52SBarry Smith Output Parameter: 3891db36a52SBarry Smith . n - number of non-redundant values 3901db36a52SBarry Smith 3911db36a52SBarry Smith Level: intermediate 3921db36a52SBarry Smith 393db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()` 3941db36a52SBarry Smith @*/ 3959371c9d4SSatish Balay PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[]) { 3961db36a52SBarry Smith PetscFunctionBegin; 3975f80ce2aSJacob Faibussowitsch PetscValidIntPointer(n, 1); 3989566063dSJacob Faibussowitsch PetscCall(PetscSortInt(*n, X)); 3999566063dSJacob Faibussowitsch PetscCall(PetscSortedRemoveDupsInt(n, X)); 4001db36a52SBarry Smith PetscFunctionReturn(0); 4011db36a52SBarry Smith } 4021db36a52SBarry Smith 40360e03357SMatthew G Knepley /*@ 404*811af0c4SBarry Smith PetscFindInt - Finds `PetscInt` in a sorted array of `PetscInt` 40560e03357SMatthew G Knepley 40660e03357SMatthew G Knepley Not Collective 40760e03357SMatthew G Knepley 40860e03357SMatthew G Knepley Input Parameters: 40960e03357SMatthew G Knepley + key - the integer to locate 410d9f1162dSJed Brown . n - number of values in the array 41159e16d97SJunchao Zhang - X - array of integers 41260e03357SMatthew G Knepley 41360e03357SMatthew G Knepley Output Parameter: 414d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 41560e03357SMatthew G Knepley 41660e03357SMatthew G Knepley Level: intermediate 41760e03357SMatthew G Knepley 418db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()` 41960e03357SMatthew G Knepley @*/ 4209371c9d4SSatish Balay PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc) { 421d9f1162dSJed Brown PetscInt lo = 0, hi = n; 42260e03357SMatthew G Knepley 42360e03357SMatthew G Knepley PetscFunctionBegin; 4245f80ce2aSJacob Faibussowitsch PetscValidIntPointer(loc, 4); 4259371c9d4SSatish Balay if (!n) { 4269371c9d4SSatish Balay *loc = -1; 4279371c9d4SSatish Balay PetscFunctionReturn(0); 4289371c9d4SSatish Balay } 429dadcf809SJacob Faibussowitsch PetscValidIntPointer(X, 3); 4306a7c706bSVaclav Hapla PetscCheckSorted(n, X); 431d9f1162dSJed Brown while (hi - lo > 1) { 432d9f1162dSJed Brown PetscInt mid = lo + (hi - lo) / 2; 43359e16d97SJunchao Zhang if (key < X[mid]) hi = mid; 434d9f1162dSJed Brown else lo = mid; 43560e03357SMatthew G Knepley } 43659e16d97SJunchao Zhang *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 43760e03357SMatthew G Knepley PetscFunctionReturn(0); 43860e03357SMatthew G Knepley } 43960e03357SMatthew G Knepley 440d2aeb606SJed Brown /*@ 441*811af0c4SBarry Smith PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates 442f1cab4e1SJunchao Zhang 443f1cab4e1SJunchao Zhang Not Collective 444f1cab4e1SJunchao Zhang 445f1cab4e1SJunchao Zhang Input Parameters: 446f1cab4e1SJunchao Zhang + n - number of values in the array 447f1cab4e1SJunchao Zhang - X - array of integers 448f1cab4e1SJunchao Zhang 449f1cab4e1SJunchao Zhang Output Parameter: 450f1cab4e1SJunchao Zhang . dups - True if the array has dups, otherwise false 451f1cab4e1SJunchao Zhang 452f1cab4e1SJunchao Zhang Level: intermediate 453f1cab4e1SJunchao Zhang 454db781477SPatrick Sanan .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()` 455f1cab4e1SJunchao Zhang @*/ 4569371c9d4SSatish Balay PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups) { 457f1cab4e1SJunchao Zhang PetscInt i; 458f1cab4e1SJunchao Zhang PetscHSetI ht; 459f1cab4e1SJunchao Zhang PetscBool missing; 460f1cab4e1SJunchao Zhang 461f1cab4e1SJunchao Zhang PetscFunctionBegin; 4625f80ce2aSJacob Faibussowitsch if (n) PetscValidIntPointer(X, 2); 4635f80ce2aSJacob Faibussowitsch PetscValidBoolPointer(dups, 3); 464f1cab4e1SJunchao Zhang *dups = PETSC_FALSE; 465f1cab4e1SJunchao Zhang if (n > 1) { 4669566063dSJacob Faibussowitsch PetscCall(PetscHSetICreate(&ht)); 4679566063dSJacob Faibussowitsch PetscCall(PetscHSetIResize(ht, n)); 468f1cab4e1SJunchao Zhang for (i = 0; i < n; i++) { 4699566063dSJacob Faibussowitsch PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing)); 4709371c9d4SSatish Balay if (!missing) { 4719371c9d4SSatish Balay *dups = PETSC_TRUE; 4729371c9d4SSatish Balay break; 4739371c9d4SSatish Balay } 474f1cab4e1SJunchao Zhang } 4759566063dSJacob Faibussowitsch PetscCall(PetscHSetIDestroy(&ht)); 476f1cab4e1SJunchao Zhang } 477f1cab4e1SJunchao Zhang PetscFunctionReturn(0); 478f1cab4e1SJunchao Zhang } 479f1cab4e1SJunchao Zhang 480f1cab4e1SJunchao Zhang /*@ 481*811af0c4SBarry Smith PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt` 482d2aeb606SJed Brown 483d2aeb606SJed Brown Not Collective 484d2aeb606SJed Brown 485d2aeb606SJed Brown Input Parameters: 486d2aeb606SJed Brown + key - the integer to locate 487d2aeb606SJed Brown . n - number of values in the array 48859e16d97SJunchao Zhang - X - array of integers 489d2aeb606SJed Brown 490d2aeb606SJed Brown Output Parameter: 491d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 492d2aeb606SJed Brown 493d2aeb606SJed Brown Level: intermediate 494d2aeb606SJed Brown 495db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()` 496d2aeb606SJed Brown @*/ 4979371c9d4SSatish Balay PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc) { 498d2aeb606SJed Brown PetscInt lo = 0, hi = n; 499d2aeb606SJed Brown 500d2aeb606SJed Brown PetscFunctionBegin; 501dadcf809SJacob Faibussowitsch PetscValidIntPointer(loc, 4); 5029371c9d4SSatish Balay if (!n) { 5039371c9d4SSatish Balay *loc = -1; 5049371c9d4SSatish Balay PetscFunctionReturn(0); 5059371c9d4SSatish Balay } 506dadcf809SJacob Faibussowitsch PetscValidIntPointer(X, 3); 5076a7c706bSVaclav Hapla PetscCheckSorted(n, X); 508d2aeb606SJed Brown while (hi - lo > 1) { 509d2aeb606SJed Brown PetscInt mid = lo + (hi - lo) / 2; 51059e16d97SJunchao Zhang if (key < X[mid]) hi = mid; 511d2aeb606SJed Brown else lo = mid; 512d2aeb606SJed Brown } 51359e16d97SJunchao Zhang *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 514e5c89e4eSSatish Balay PetscFunctionReturn(0); 515e5c89e4eSSatish Balay } 516e5c89e4eSSatish Balay 517e5c89e4eSSatish Balay /*@ 518*811af0c4SBarry Smith PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order; 519*811af0c4SBarry Smith changes a second array of `PetscInt` to match the sorted first array. 520e5c89e4eSSatish Balay 521e5c89e4eSSatish Balay Not Collective 522e5c89e4eSSatish Balay 523e5c89e4eSSatish Balay Input Parameters: 524e5c89e4eSSatish Balay + n - number of values 52559e16d97SJunchao Zhang . X - array of integers 52659e16d97SJunchao Zhang - Y - second array of integers 527e5c89e4eSSatish Balay 528e5c89e4eSSatish Balay Level: intermediate 529e5c89e4eSSatish Balay 530db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()` 531e5c89e4eSSatish Balay @*/ 5329371c9d4SSatish Balay PetscErrorCode PetscSortIntWithArray(PetscInt n, PetscInt X[], PetscInt Y[]) { 53359e16d97SJunchao Zhang PetscInt pivot, t1, t2; 534e5c89e4eSSatish Balay 535e5c89e4eSSatish Balay PetscFunctionBegin; 5365f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2); 537c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 538c1f0200aSDmitry Karpeev } 539c1f0200aSDmitry Karpeev 540c1f0200aSDmitry Karpeev /*@ 541*811af0c4SBarry Smith PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order; 542*811af0c4SBarry Smith changes a pair of `PetscInt` arrays to match the sorted first array. 543c1f0200aSDmitry Karpeev 544c1f0200aSDmitry Karpeev Not Collective 545c1f0200aSDmitry Karpeev 546c1f0200aSDmitry Karpeev Input Parameters: 547c1f0200aSDmitry Karpeev + n - number of values 54859e16d97SJunchao Zhang . X - array of integers 54959e16d97SJunchao Zhang . Y - second array of integers (first array of the pair) 55059e16d97SJunchao Zhang - Z - third array of integers (second array of the pair) 551c1f0200aSDmitry Karpeev 552c1f0200aSDmitry Karpeev Level: intermediate 553c1f0200aSDmitry Karpeev 554db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()` 555c1f0200aSDmitry Karpeev @*/ 5569371c9d4SSatish Balay PetscErrorCode PetscSortIntWithArrayPair(PetscInt n, PetscInt X[], PetscInt Y[], PetscInt Z[]) { 55759e16d97SJunchao Zhang PetscInt pivot, t1, t2, t3; 558c1f0200aSDmitry Karpeev 559c1f0200aSDmitry Karpeev PetscFunctionBegin; 5605f80ce2aSJacob Faibussowitsch QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3); 561c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 562c1f0200aSDmitry Karpeev } 563c1f0200aSDmitry Karpeev 56417d7d925SStefano Zampini /*@ 565*811af0c4SBarry Smith PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order; 566*811af0c4SBarry Smith changes a second array of `PetscCount` to match the sorted first array. 567981bb840SJunchao Zhang 568981bb840SJunchao Zhang Not Collective 569981bb840SJunchao Zhang 570981bb840SJunchao Zhang Input Parameters: 571981bb840SJunchao Zhang + n - number of values 572981bb840SJunchao Zhang . X - array of integers 573981bb840SJunchao Zhang - Y - second array of PetscCounts (signed integers) 574981bb840SJunchao Zhang 575981bb840SJunchao Zhang Level: intermediate 576981bb840SJunchao Zhang 577db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 578981bb840SJunchao Zhang @*/ 5799371c9d4SSatish Balay PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[]) { 580981bb840SJunchao Zhang PetscInt pivot, t1; 581981bb840SJunchao Zhang PetscCount t2; 582981bb840SJunchao Zhang 583981bb840SJunchao Zhang PetscFunctionBegin; 5845f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2); 585981bb840SJunchao Zhang PetscFunctionReturn(0); 586981bb840SJunchao Zhang } 587981bb840SJunchao Zhang 588981bb840SJunchao Zhang /*@ 589*811af0c4SBarry Smith PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order; 590*811af0c4SBarry Smith changes a `PetscInt` array and a `PetscCount` array to match the sorted first array. 591981bb840SJunchao Zhang 592981bb840SJunchao Zhang Not Collective 593981bb840SJunchao Zhang 594981bb840SJunchao Zhang Input Parameters: 595981bb840SJunchao Zhang + n - number of values 596981bb840SJunchao Zhang . X - array of integers 597981bb840SJunchao Zhang . Y - second array of integers (first array of the pair) 598981bb840SJunchao Zhang - Z - third array of PetscCounts (second array of the pair) 599981bb840SJunchao Zhang 600981bb840SJunchao Zhang Level: intermediate 601981bb840SJunchao Zhang 602*811af0c4SBarry Smith Note: 603981bb840SJunchao 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. 604981bb840SJunchao Zhang 605db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()` 606981bb840SJunchao Zhang @*/ 6079371c9d4SSatish Balay PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[]) { 608981bb840SJunchao Zhang PetscInt pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */ 609981bb840SJunchao Zhang PetscCount t3; /* temp for Z[] */ 610981bb840SJunchao Zhang 611981bb840SJunchao Zhang PetscFunctionBegin; 6125f80ce2aSJacob Faibussowitsch QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3); 613981bb840SJunchao Zhang PetscFunctionReturn(0); 614981bb840SJunchao Zhang } 615981bb840SJunchao Zhang 616981bb840SJunchao Zhang /*@ 617*811af0c4SBarry Smith PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted. 6186a7c706bSVaclav Hapla 6196a7c706bSVaclav Hapla Not Collective 6206a7c706bSVaclav Hapla 6216a7c706bSVaclav Hapla Input Parameters: 6226a7c706bSVaclav Hapla + n - number of values 6236a7c706bSVaclav Hapla - X - array of integers 6246a7c706bSVaclav Hapla 6256a7c706bSVaclav Hapla Output Parameters: 6266a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 6276a7c706bSVaclav Hapla 6286a7c706bSVaclav Hapla Level: intermediate 6296a7c706bSVaclav Hapla 630db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()` 6316a7c706bSVaclav Hapla @*/ 6329371c9d4SSatish Balay PetscErrorCode PetscSortedMPIInt(PetscInt n, const PetscMPIInt X[], PetscBool *sorted) { 6336a7c706bSVaclav Hapla PetscFunctionBegin; 6346a7c706bSVaclav Hapla PetscSorted(n, X, *sorted); 6356a7c706bSVaclav Hapla PetscFunctionReturn(0); 6366a7c706bSVaclav Hapla } 6376a7c706bSVaclav Hapla 6386a7c706bSVaclav Hapla /*@ 639*811af0c4SBarry Smith PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order. 64017d7d925SStefano Zampini 64117d7d925SStefano Zampini Not Collective 64217d7d925SStefano Zampini 64317d7d925SStefano Zampini Input Parameters: 64417d7d925SStefano Zampini + n - number of values 64559e16d97SJunchao Zhang - X - array of integers 64617d7d925SStefano Zampini 64717d7d925SStefano Zampini Level: intermediate 64817d7d925SStefano Zampini 649*811af0c4SBarry Smith Note: 650676f2a66SJacob Faibussowitsch This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array 651a5b23f4aSJose E. Roman is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their 652676f2a66SJacob Faibussowitsch code to see which routine is fastest. 653676f2a66SJacob Faibussowitsch 654db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()` 65517d7d925SStefano Zampini @*/ 6569371c9d4SSatish Balay PetscErrorCode PetscSortMPIInt(PetscInt n, PetscMPIInt X[]) { 65759e16d97SJunchao Zhang PetscMPIInt pivot, t1; 65817d7d925SStefano Zampini 65917d7d925SStefano Zampini PetscFunctionBegin; 6605f80ce2aSJacob Faibussowitsch QuickSort1(PetscSortMPIInt, X, n, pivot, t1); 66117d7d925SStefano Zampini PetscFunctionReturn(0); 66217d7d925SStefano Zampini } 66317d7d925SStefano Zampini 66417d7d925SStefano Zampini /*@ 665*811af0c4SBarry Smith PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries 66617d7d925SStefano Zampini 66717d7d925SStefano Zampini Not Collective 66817d7d925SStefano Zampini 66917d7d925SStefano Zampini Input Parameters: 67017d7d925SStefano Zampini + n - number of values 67159e16d97SJunchao Zhang - X - array of integers 67217d7d925SStefano Zampini 67317d7d925SStefano Zampini Output Parameter: 67417d7d925SStefano Zampini . n - number of non-redundant values 67517d7d925SStefano Zampini 67617d7d925SStefano Zampini Level: intermediate 67717d7d925SStefano Zampini 678db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()` 67917d7d925SStefano Zampini @*/ 6809371c9d4SSatish Balay PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[]) { 6815f80ce2aSJacob Faibussowitsch PetscInt s = 0, N = *n, b = 0; 68217d7d925SStefano Zampini 68317d7d925SStefano Zampini PetscFunctionBegin; 6849566063dSJacob Faibussowitsch PetscCall(PetscSortMPIInt(N, X)); 6855f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < N - 1; i++) { 68659e16d97SJunchao Zhang if (X[b + s + 1] != X[b]) { 6879371c9d4SSatish Balay X[b + 1] = X[b + s + 1]; 6889371c9d4SSatish Balay b++; 689a297a907SKarl Rupp } else s++; 69017d7d925SStefano Zampini } 69117d7d925SStefano Zampini *n = N - s; 69217d7d925SStefano Zampini PetscFunctionReturn(0); 69317d7d925SStefano Zampini } 694c1f0200aSDmitry Karpeev 6954d615ea0SBarry Smith /*@ 696*811af0c4SBarry Smith PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order; 697*811af0c4SBarry Smith changes a second `PetscMPIInt` array to match the sorted first array. 6984d615ea0SBarry Smith 6994d615ea0SBarry Smith Not Collective 7004d615ea0SBarry Smith 7014d615ea0SBarry Smith Input Parameters: 7024d615ea0SBarry Smith + n - number of values 70359e16d97SJunchao Zhang . X - array of integers 70459e16d97SJunchao Zhang - Y - second array of integers 7054d615ea0SBarry Smith 7064d615ea0SBarry Smith Level: intermediate 7074d615ea0SBarry Smith 708db781477SPatrick Sanan .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()` 7094d615ea0SBarry Smith @*/ 7109371c9d4SSatish Balay PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n, PetscMPIInt X[], PetscMPIInt Y[]) { 71159e16d97SJunchao Zhang PetscMPIInt pivot, t1, t2; 7124d615ea0SBarry Smith 7134d615ea0SBarry Smith PetscFunctionBegin; 7145f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2); 7155569a785SJunchao Zhang PetscFunctionReturn(0); 7165569a785SJunchao Zhang } 7175569a785SJunchao Zhang 7185569a785SJunchao Zhang /*@ 719*811af0c4SBarry Smith PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order; 720*811af0c4SBarry Smith changes a second array of `PetscInt` to match the sorted first array. 7215569a785SJunchao Zhang 7225569a785SJunchao Zhang Not Collective 7235569a785SJunchao Zhang 7245569a785SJunchao Zhang Input Parameters: 7255569a785SJunchao Zhang + n - number of values 72659e16d97SJunchao Zhang . X - array of MPI integers 72759e16d97SJunchao Zhang - Y - second array of Petsc integers 7285569a785SJunchao Zhang 7295569a785SJunchao Zhang Level: intermediate 7305569a785SJunchao Zhang 731*811af0c4SBarry Smith Note: 732*811af0c4SBarry Smith This routine is useful when one needs to sort MPI ranks with other integer arrays. 7335569a785SJunchao Zhang 734db781477SPatrick Sanan .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()` 7355569a785SJunchao Zhang @*/ 7369371c9d4SSatish Balay PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n, PetscMPIInt X[], PetscInt Y[]) { 73759e16d97SJunchao Zhang PetscMPIInt pivot, t1; 7385569a785SJunchao Zhang PetscInt t2; 7395569a785SJunchao Zhang 7405569a785SJunchao Zhang PetscFunctionBegin; 7415f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2); 742e5c89e4eSSatish Balay PetscFunctionReturn(0); 743e5c89e4eSSatish Balay } 744e5c89e4eSSatish Balay 745e5c89e4eSSatish Balay /*@ 746*811af0c4SBarry Smith PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order; 747*811af0c4SBarry Smith changes a second `PetscScalar` array to match the sorted first array. 748e5c89e4eSSatish Balay 749e5c89e4eSSatish Balay Not Collective 750e5c89e4eSSatish Balay 751e5c89e4eSSatish Balay Input Parameters: 752e5c89e4eSSatish Balay + n - number of values 75359e16d97SJunchao Zhang . X - array of integers 75459e16d97SJunchao Zhang - Y - second array of scalars 755e5c89e4eSSatish Balay 756e5c89e4eSSatish Balay Level: intermediate 757e5c89e4eSSatish Balay 758db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 759e5c89e4eSSatish Balay @*/ 7609371c9d4SSatish Balay PetscErrorCode PetscSortIntWithScalarArray(PetscInt n, PetscInt X[], PetscScalar Y[]) { 76159e16d97SJunchao Zhang PetscInt pivot, t1; 76259e16d97SJunchao Zhang PetscScalar t2; 763e5c89e4eSSatish Balay 764e5c89e4eSSatish Balay PetscFunctionBegin; 7655f80ce2aSJacob Faibussowitsch QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2); 76617df18f2SToby Isaac PetscFunctionReturn(0); 76717df18f2SToby Isaac } 76817df18f2SToby Isaac 7696c2863d0SBarry Smith /*@C 770*811af0c4SBarry Smith PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order; 77117df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must 77217df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting. 77317df18f2SToby Isaac 77417df18f2SToby Isaac Not Collective 77517df18f2SToby Isaac 77617df18f2SToby Isaac Input Parameters: 77717df18f2SToby Isaac + n - number of values 77859e16d97SJunchao Zhang . X - array of integers 77959e16d97SJunchao Zhang . Y - second array of data 78017df18f2SToby Isaac . size - sizeof elements in the data array in bytes 78159e16d97SJunchao Zhang - t2 - workspace of "size" bytes used when sorting 78217df18f2SToby Isaac 78317df18f2SToby Isaac Level: intermediate 78417df18f2SToby Isaac 785db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 78617df18f2SToby Isaac @*/ 7879371c9d4SSatish Balay PetscErrorCode PetscSortIntWithDataArray(PetscInt n, PetscInt X[], void *Y, size_t size, void *t2) { 78859e16d97SJunchao Zhang char *YY = (char *)Y; 7895f80ce2aSJacob Faibussowitsch PetscInt t1, pivot, hi = n - 1; 79017df18f2SToby Isaac 79117df18f2SToby Isaac PetscFunctionBegin; 79217df18f2SToby Isaac if (n < 8) { 7935f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < n; i++) { 79459e16d97SJunchao Zhang pivot = X[i]; 7955f80ce2aSJacob Faibussowitsch for (PetscInt j = i + 1; j < n; j++) { 79659e16d97SJunchao Zhang if (pivot > X[j]) { 79759e16d97SJunchao Zhang SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size); 79859e16d97SJunchao Zhang pivot = X[i]; 79917df18f2SToby Isaac } 80017df18f2SToby Isaac } 80117df18f2SToby Isaac } 80217df18f2SToby Isaac } else { 80359e16d97SJunchao Zhang /* Two way partition */ 8045f80ce2aSJacob Faibussowitsch PetscInt l = 0, r = hi; 8055f80ce2aSJacob Faibussowitsch 8065f80ce2aSJacob Faibussowitsch pivot = X[MEDIAN(X, hi)]; 80759e16d97SJunchao Zhang while (1) { 80859e16d97SJunchao Zhang while (X[l] < pivot) l++; 80959e16d97SJunchao Zhang while (X[r] > pivot) r--; 8109371c9d4SSatish Balay if (l >= r) { 8119371c9d4SSatish Balay r++; 8129371c9d4SSatish Balay break; 8139371c9d4SSatish Balay } 81459e16d97SJunchao Zhang SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size); 81559e16d97SJunchao Zhang l++; 81659e16d97SJunchao Zhang r--; 81759e16d97SJunchao Zhang } 8189566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2)); 8199566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2)); 82017df18f2SToby Isaac } 82117df18f2SToby Isaac PetscFunctionReturn(0); 82217df18f2SToby Isaac } 82317df18f2SToby Isaac 82421e72a00SBarry Smith /*@ 825*811af0c4SBarry Smith PetscMergeIntArray - Merges two SORTED `PetscInt` arrays, removes duplicate elements. 82621e72a00SBarry Smith 82721e72a00SBarry Smith Not Collective 82821e72a00SBarry Smith 82921e72a00SBarry Smith Input Parameters: 83021e72a00SBarry Smith + an - number of values in the first array 83121e72a00SBarry Smith . aI - first sorted array of integers 83221e72a00SBarry Smith . bn - number of values in the second array 83321e72a00SBarry Smith - bI - second array of integers 83421e72a00SBarry Smith 83521e72a00SBarry Smith Output Parameters: 83621e72a00SBarry Smith + n - number of values in the merged array 8376c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided 83821e72a00SBarry Smith 83921e72a00SBarry Smith Level: intermediate 84021e72a00SBarry Smith 841db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 84221e72a00SBarry Smith @*/ 8439371c9d4SSatish Balay PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) { 84421e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k; 84521e72a00SBarry Smith 846362febeeSStefano Zampini PetscFunctionBegin; 84721e72a00SBarry Smith if (!L_) { 8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(an + bn, L)); 84921e72a00SBarry Smith L_ = *L; 85021e72a00SBarry Smith } 85121e72a00SBarry Smith k = ak = bk = 0; 85221e72a00SBarry Smith while (ak < an && bk < bn) { 85321e72a00SBarry Smith if (aI[ak] == bI[bk]) { 85421e72a00SBarry Smith L_[k] = aI[ak]; 85521e72a00SBarry Smith ++ak; 85621e72a00SBarry Smith ++bk; 85721e72a00SBarry Smith ++k; 85821e72a00SBarry Smith } else if (aI[ak] < bI[bk]) { 85921e72a00SBarry Smith L_[k] = aI[ak]; 86021e72a00SBarry Smith ++ak; 86121e72a00SBarry Smith ++k; 86221e72a00SBarry Smith } else { 86321e72a00SBarry Smith L_[k] = bI[bk]; 86421e72a00SBarry Smith ++bk; 86521e72a00SBarry Smith ++k; 86621e72a00SBarry Smith } 86721e72a00SBarry Smith } 86821e72a00SBarry Smith if (ak < an) { 8699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak)); 87021e72a00SBarry Smith k += (an - ak); 87121e72a00SBarry Smith } 87221e72a00SBarry Smith if (bk < bn) { 8739566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk)); 87421e72a00SBarry Smith k += (bn - bk); 87521e72a00SBarry Smith } 87621e72a00SBarry Smith *n = k; 87721e72a00SBarry Smith PetscFunctionReturn(0); 87821e72a00SBarry Smith } 879b4301105SBarry Smith 880c1f0200aSDmitry Karpeev /*@ 881*811af0c4SBarry Smith PetscMergeIntArrayPair - Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`. 882c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 883c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 884c1f0200aSDmitry Karpeev 885c1f0200aSDmitry Karpeev Not Collective 886c1f0200aSDmitry Karpeev 887c1f0200aSDmitry Karpeev Input Parameters: 888c1f0200aSDmitry Karpeev + an - number of values in the first array 889c1f0200aSDmitry Karpeev . aI - first sorted array of integers 890c1f0200aSDmitry Karpeev . aJ - first additional array of integers 891c1f0200aSDmitry Karpeev . bn - number of values in the second array 892c1f0200aSDmitry Karpeev . bI - second array of integers 893c1f0200aSDmitry Karpeev - bJ - second additional of integers 894c1f0200aSDmitry Karpeev 895c1f0200aSDmitry Karpeev Output Parameters: 896c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 89714c49006SJed Brown . L - merged sorted array 898c1f0200aSDmitry Karpeev - J - merged additional array 899c1f0200aSDmitry Karpeev 900*811af0c4SBarry Smith Note: 90195452b02SPatrick 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 902*811af0c4SBarry Smith 903c1f0200aSDmitry Karpeev Level: intermediate 904c1f0200aSDmitry Karpeev 905db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 906c1f0200aSDmitry Karpeev @*/ 9079371c9d4SSatish Balay PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J) { 90870d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k; 909c1f0200aSDmitry Karpeev 91014c49006SJed Brown PetscFunctionBegin; 911dadcf809SJacob Faibussowitsch PetscValidPointer(L, 8); 912dadcf809SJacob Faibussowitsch PetscValidPointer(J, 9); 913c1f0200aSDmitry Karpeev n_ = an + bn; 914c1f0200aSDmitry Karpeev *n = n_; 91548a46eb9SPierre Jolivet if (!*L) PetscCall(PetscMalloc1(n_, L)); 916d7aa01a8SHong Zhang L_ = *L; 91748a46eb9SPierre Jolivet if (!*J) PetscCall(PetscMalloc1(n_, J)); 918c1f0200aSDmitry Karpeev J_ = *J; 919c1f0200aSDmitry Karpeev k = ak = bk = 0; 920c1f0200aSDmitry Karpeev while (ak < an && bk < bn) { 921c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 922d7aa01a8SHong Zhang L_[k] = aI[ak]; 923c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 924c1f0200aSDmitry Karpeev ++ak; 925c1f0200aSDmitry Karpeev ++k; 926a297a907SKarl Rupp } else { 927d7aa01a8SHong Zhang L_[k] = bI[bk]; 928c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 929c1f0200aSDmitry Karpeev ++bk; 930c1f0200aSDmitry Karpeev ++k; 931c1f0200aSDmitry Karpeev } 932c1f0200aSDmitry Karpeev } 933c1f0200aSDmitry Karpeev if (ak < an) { 9349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak)); 9359566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak)); 936c1f0200aSDmitry Karpeev k += (an - ak); 937c1f0200aSDmitry Karpeev } 938c1f0200aSDmitry Karpeev if (bk < bn) { 9399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk)); 9409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk)); 941c1f0200aSDmitry Karpeev } 942c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 943c1f0200aSDmitry Karpeev } 944c1f0200aSDmitry Karpeev 945e498c390SJed Brown /*@ 946*811af0c4SBarry Smith PetscMergeMPIIntArray - Merges two SORTED `PetscMPIInt` arrays. 947e498c390SJed Brown 948e498c390SJed Brown Not Collective 949e498c390SJed Brown 950e498c390SJed Brown Input Parameters: 951e498c390SJed Brown + an - number of values in the first array 952e498c390SJed Brown . aI - first sorted array of integers 953e498c390SJed Brown . bn - number of values in the second array 954e498c390SJed Brown - bI - second array of integers 955e498c390SJed Brown 956e498c390SJed Brown Output Parameters: 957e498c390SJed Brown + n - number of values in the merged array (<= an + bn) 958e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed 959e498c390SJed Brown 960e498c390SJed Brown Level: intermediate 961e498c390SJed Brown 962db781477SPatrick Sanan .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()` 963e498c390SJed Brown @*/ 9649371c9d4SSatish Balay PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L) { 965e498c390SJed Brown PetscInt ai, bi, k; 966e498c390SJed Brown 967e498c390SJed Brown PetscFunctionBegin; 9689566063dSJacob Faibussowitsch if (!*L) PetscCall(PetscMalloc1((an + bn), L)); 969e498c390SJed Brown for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) { 970e498c390SJed Brown PetscInt t = -1; 971e498c390SJed Brown for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai]; 9729371c9d4SSatish Balay for (; bi < bn && bI[bi] == t; bi++) 9739371c9d4SSatish Balay ; 974e498c390SJed Brown for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi]; 9759371c9d4SSatish Balay for (; ai < an && aI[ai] == t; ai++) 9769371c9d4SSatish Balay ; 977e498c390SJed Brown } 978e498c390SJed Brown *n = k; 979e498c390SJed Brown PetscFunctionReturn(0); 980e498c390SJed Brown } 981e5c89e4eSSatish Balay 9826c2863d0SBarry Smith /*@C 98342eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 98442eaa7f3SBarry Smith 98542eaa7f3SBarry Smith Not Collective 98642eaa7f3SBarry Smith 98742eaa7f3SBarry Smith Input Parameters: 98842eaa7f3SBarry Smith + n - number of values 98942eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 99042eaa7f3SBarry Smith - parentid - indicates the parent of each entry 99142eaa7f3SBarry Smith 99242eaa7f3SBarry Smith Output Parameters: 99342eaa7f3SBarry Smith + Nlevels - the number of levels 99442eaa7f3SBarry Smith . Level - for each node tells its level 99542eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 99642eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 99742eaa7f3SBarry Smith - Column - for each id tells its column index 99842eaa7f3SBarry Smith 9996c2863d0SBarry Smith Level: developer 100042eaa7f3SBarry Smith 1001*811af0c4SBarry Smith Note: 100295452b02SPatrick Sanan This code is not currently used 100342eaa7f3SBarry Smith 1004db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()` 100542eaa7f3SBarry Smith @*/ 10069371c9d4SSatish Balay PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column) { 100742eaa7f3SBarry Smith PetscInt i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column; 1008ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 100942eaa7f3SBarry Smith 101042eaa7f3SBarry Smith PetscFunctionBegin; 101108401ef6SPierre Jolivet PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set"); 101242eaa7f3SBarry Smith for (i = 0; i < n; i++) { 101342eaa7f3SBarry Smith if (mask[i]) continue; 101408401ef6SPierre Jolivet PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent"); 1015cc73adaaSBarry Smith PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked"); 101642eaa7f3SBarry Smith } 101742eaa7f3SBarry Smith 101842eaa7f3SBarry Smith for (i = 0; i < n; i++) { 101942eaa7f3SBarry Smith if (!mask[i]) nmask++; 102042eaa7f3SBarry Smith } 102142eaa7f3SBarry Smith 102242eaa7f3SBarry Smith /* determine the level in the tree of each node */ 10239566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n, &level)); 1024a297a907SKarl Rupp 102542eaa7f3SBarry Smith level[0] = 1; 102642eaa7f3SBarry Smith while (!done) { 102742eaa7f3SBarry Smith done = PETSC_TRUE; 102842eaa7f3SBarry Smith for (i = 0; i < n; i++) { 102942eaa7f3SBarry Smith if (mask[i]) continue; 103042eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 103142eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 103242eaa7f3SBarry Smith } 103342eaa7f3SBarry Smith } 103442eaa7f3SBarry Smith for (i = 0; i < n; i++) { 103542eaa7f3SBarry Smith level[i]--; 103642eaa7f3SBarry Smith nlevels = PetscMax(nlevels, level[i]); 103742eaa7f3SBarry Smith } 103842eaa7f3SBarry Smith 103942eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 10409566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nlevels, &levelcnt)); 104142eaa7f3SBarry Smith for (i = 0; i < n; i++) { 104242eaa7f3SBarry Smith if (mask[i]) continue; 104342eaa7f3SBarry Smith levelcnt[level[i] - 1]++; 104442eaa7f3SBarry Smith } 1045a297a907SKarl Rupp for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]); 104642eaa7f3SBarry Smith 104742eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 10489566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid)); 10499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nmask, &idbylevel)); 105042eaa7f3SBarry Smith for (j = 1; j <= nlevels; j++) { 105142eaa7f3SBarry Smith cnt = 0; 105242eaa7f3SBarry Smith for (i = 0; i < n; i++) { 105342eaa7f3SBarry Smith if (mask[i]) continue; 105442eaa7f3SBarry Smith if (level[i] != j) continue; 105542eaa7f3SBarry Smith workid[cnt] = i; 105642eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 105742eaa7f3SBarry Smith } 105842eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 105942eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 10609566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(cnt,workparentid,workid)); 106142eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 106242eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 10639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt)); 106442eaa7f3SBarry Smith tcnt += cnt; 106542eaa7f3SBarry Smith } 106608401ef6SPierre Jolivet PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes"); 10679566063dSJacob Faibussowitsch PetscCall(PetscFree2(workid, workparentid)); 106842eaa7f3SBarry Smith 106942eaa7f3SBarry Smith /* for each node list its column */ 10709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &column)); 107142eaa7f3SBarry Smith cnt = 0; 107242eaa7f3SBarry Smith for (j = 0; j < nlevels; j++) { 1073ad540459SPierre Jolivet for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i; 107442eaa7f3SBarry Smith } 107542eaa7f3SBarry Smith 107642eaa7f3SBarry Smith *Nlevels = nlevels; 107742eaa7f3SBarry Smith *Level = level; 107842eaa7f3SBarry Smith *Levelcnt = levelcnt; 107942eaa7f3SBarry Smith *Idbylevel = idbylevel; 108042eaa7f3SBarry Smith *Column = column; 108142eaa7f3SBarry Smith PetscFunctionReturn(0); 108242eaa7f3SBarry Smith } 1083ce605777SToby Isaac 1084ce605777SToby Isaac /*@ 1085*811af0c4SBarry Smith PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted. 1086ce605777SToby Isaac 1087ce605777SToby Isaac Collective 1088ce605777SToby Isaac 1089ce605777SToby Isaac Input Parameters: 1090ce605777SToby Isaac + comm - the MPI communicator 1091ce605777SToby Isaac . n - the local number of integers 1092ce605777SToby Isaac - keys - the local array of integers 1093ce605777SToby Isaac 1094ce605777SToby Isaac Output Parameters: 1095ce605777SToby Isaac . is_sorted - whether the array is globally sorted 1096ce605777SToby Isaac 1097ce605777SToby Isaac Level: developer 1098ce605777SToby Isaac 1099db781477SPatrick Sanan .seealso: `PetscParallelSortInt()` 1100ce605777SToby Isaac @*/ 11019371c9d4SSatish Balay PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted) { 1102ce605777SToby Isaac PetscBool sorted; 1103ce605777SToby Isaac PetscInt i, min, max, prevmax; 11042a1da528SToby Isaac PetscMPIInt rank; 1105ce605777SToby Isaac 1106ce605777SToby Isaac PetscFunctionBegin; 1107ce605777SToby Isaac sorted = PETSC_TRUE; 1108ce605777SToby Isaac min = PETSC_MAX_INT; 1109ce605777SToby Isaac max = PETSC_MIN_INT; 1110ce605777SToby Isaac if (n) { 1111ce605777SToby Isaac min = keys[0]; 1112ce605777SToby Isaac max = keys[0]; 1113ce605777SToby Isaac } 1114ce605777SToby Isaac for (i = 1; i < n; i++) { 1115ce605777SToby Isaac if (keys[i] < keys[i - 1]) break; 1116ce605777SToby Isaac min = PetscMin(min, keys[i]); 1117ce605777SToby Isaac max = PetscMax(max, keys[i]); 1118ce605777SToby Isaac } 1119ce605777SToby Isaac if (i < n) sorted = PETSC_FALSE; 1120ce605777SToby Isaac prevmax = PETSC_MIN_INT; 11219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm)); 11229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1123dd400576SPatrick Sanan if (rank == 0) prevmax = PETSC_MIN_INT; 1124ce605777SToby Isaac if (prevmax > min) sorted = PETSC_FALSE; 11259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm)); 1126ce605777SToby Isaac PetscFunctionReturn(0); 1127ce605777SToby Isaac } 1128