1676f2a66SJacob Faibussowitsch #include <petscsys.h> /*I "petscsys.h" I*/ 2676f2a66SJacob Faibussowitsch #include <petsc/private/petscimpl.h> 3676f2a66SJacob Faibussowitsch 4d71ae5a4SJacob Faibussowitsch static inline int Compare_PetscMPIInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) 5d71ae5a4SJacob Faibussowitsch { 6676f2a66SJacob Faibussowitsch PetscMPIInt l = *(PetscMPIInt *)left, r = *(PetscMPIInt *)right; 7676f2a66SJacob Faibussowitsch return (l < r) ? -1 : (l > r); 8676f2a66SJacob Faibussowitsch } 9676f2a66SJacob Faibussowitsch 10d71ae5a4SJacob Faibussowitsch static inline int Compare_PetscInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) 11d71ae5a4SJacob Faibussowitsch { 12676f2a66SJacob Faibussowitsch PetscInt l = *(PetscInt *)left, r = *(PetscInt *)right; 13676f2a66SJacob Faibussowitsch return (l < r) ? -1 : (l > r); 14676f2a66SJacob Faibussowitsch } 15676f2a66SJacob Faibussowitsch 16d71ae5a4SJacob Faibussowitsch static inline int Compare_PetscReal_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) 17d71ae5a4SJacob Faibussowitsch { 18676f2a66SJacob Faibussowitsch PetscReal l = *(PetscReal *)left, r = *(PetscReal *)right; 19676f2a66SJacob Faibussowitsch return (l < r) ? -1 : (l > r); 20676f2a66SJacob Faibussowitsch } 21676f2a66SJacob Faibussowitsch 22676f2a66SJacob Faibussowitsch #define MIN_GALLOP_CONST_GLOBAL 8 23676f2a66SJacob Faibussowitsch static PetscInt MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL; 244d3610e3SJacob Faibussowitsch typedef int (*CompFunc)(const void *, const void *, void *); 25676f2a66SJacob Faibussowitsch 26676f2a66SJacob Faibussowitsch /* Mostly to force clang uses the builtins, but can't hurt to have gcc compilers also do the same */ 27676f2a66SJacob Faibussowitsch #if !defined(PETSC_USE_DEBUG) && defined(__GNUC__) 28d71ae5a4SJacob Faibussowitsch static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size) 29d71ae5a4SJacob Faibussowitsch { 30676f2a66SJacob Faibussowitsch __builtin_memcpy(t, b, size); 31676f2a66SJacob Faibussowitsch __builtin_memmove(b, a, size); 32676f2a66SJacob Faibussowitsch __builtin_memcpy(a, t, size); 33676f2a66SJacob Faibussowitsch return; 34676f2a66SJacob Faibussowitsch } 35676f2a66SJacob Faibussowitsch 36d71ae5a4SJacob Faibussowitsch static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t) 37d71ae5a4SJacob Faibussowitsch { 386e260675SJacob Faibussowitsch __builtin_memcpy(t, ar, asize); 396e260675SJacob Faibussowitsch __builtin_memmove(ar, al, asize); 406e260675SJacob Faibussowitsch __builtin_memcpy(al, t, asize); 416e260675SJacob Faibussowitsch __builtin_memcpy(t, br, bsize); 426e260675SJacob Faibussowitsch __builtin_memmove(br, bl, bsize); 436e260675SJacob Faibussowitsch __builtin_memcpy(bl, t, bsize); 44676f2a66SJacob Faibussowitsch return; 45676f2a66SJacob Faibussowitsch } 46676f2a66SJacob Faibussowitsch 47d71ae5a4SJacob Faibussowitsch static inline void Petsc_memcpy(char *dest, const char *src, size_t size) 48d71ae5a4SJacob Faibussowitsch { 49676f2a66SJacob Faibussowitsch __builtin_memcpy(dest, src, size); 50676f2a66SJacob Faibussowitsch return; 51676f2a66SJacob Faibussowitsch } 52676f2a66SJacob Faibussowitsch 53d71ae5a4SJacob Faibussowitsch static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) 54d71ae5a4SJacob Faibussowitsch { 55676f2a66SJacob Faibussowitsch __builtin_memcpy(adest, asrc, asize); 566e260675SJacob Faibussowitsch __builtin_memcpy(bdest, bsrc, bsize); 57676f2a66SJacob Faibussowitsch return; 58676f2a66SJacob Faibussowitsch } 59676f2a66SJacob Faibussowitsch 60d71ae5a4SJacob Faibussowitsch static inline void Petsc_memmove(char *dest, const char *src, size_t size) 61d71ae5a4SJacob Faibussowitsch { 62676f2a66SJacob Faibussowitsch __builtin_memmove(dest, src, size); 63676f2a66SJacob Faibussowitsch return; 64676f2a66SJacob Faibussowitsch } 65676f2a66SJacob Faibussowitsch 66d71ae5a4SJacob Faibussowitsch static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) 67d71ae5a4SJacob Faibussowitsch { 68676f2a66SJacob Faibussowitsch __builtin_memmove(adest, asrc, asize); 69676f2a66SJacob Faibussowitsch __builtin_memmove(bdest, bsrc, bsize); 70676f2a66SJacob Faibussowitsch return; 71676f2a66SJacob Faibussowitsch } 72676f2a66SJacob Faibussowitsch #else 73d71ae5a4SJacob Faibussowitsch static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size) 74d71ae5a4SJacob Faibussowitsch { 75676f2a66SJacob Faibussowitsch PetscFunctionBegin; 769566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, b, size)); 779566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(b, a, size)); 789566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(a, t, size)); 79676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 80676f2a66SJacob Faibussowitsch } 81676f2a66SJacob Faibussowitsch 82d71ae5a4SJacob Faibussowitsch static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t) 83d71ae5a4SJacob Faibussowitsch { 84676f2a66SJacob Faibussowitsch PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, ar, asize)); 869566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(ar, al, asize)); 879566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(al, t, asize)); 889566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, br, bsize)); 899566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(br, bl, bsize)); 909566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(bl, t, bsize)); 91676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 92676f2a66SJacob Faibussowitsch } 93676f2a66SJacob Faibussowitsch 94d71ae5a4SJacob Faibussowitsch static inline void Petsc_memcpy(char *dest, const char *src, size_t size) 95d71ae5a4SJacob Faibussowitsch { 96676f2a66SJacob Faibussowitsch PetscFunctionBegin; 979566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(dest, src, size)); 98676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 99676f2a66SJacob Faibussowitsch } 100676f2a66SJacob Faibussowitsch 101d71ae5a4SJacob Faibussowitsch static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) 102d71ae5a4SJacob Faibussowitsch { 103676f2a66SJacob Faibussowitsch PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(adest, asrc, asize)); 1059566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(bdest, bsrc, bsize)); 106676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 107676f2a66SJacob Faibussowitsch } 108676f2a66SJacob Faibussowitsch 109d71ae5a4SJacob Faibussowitsch static inline void Petsc_memmove(char *dest, const char *src, size_t size) 110d71ae5a4SJacob Faibussowitsch { 111676f2a66SJacob Faibussowitsch PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(dest, src, size)); 113676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 114676f2a66SJacob Faibussowitsch } 115676f2a66SJacob Faibussowitsch 116d71ae5a4SJacob Faibussowitsch static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) 117d71ae5a4SJacob Faibussowitsch { 118676f2a66SJacob Faibussowitsch PetscFunctionBegin; 1199566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(adest, asrc, asize)); 1209566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(bdest, bsrc, bsize)); 121676f2a66SJacob Faibussowitsch PetscFunctionReturnVoid(); 122676f2a66SJacob Faibussowitsch } 123676f2a66SJacob Faibussowitsch #endif 124676f2a66SJacob Faibussowitsch 125676f2a66SJacob Faibussowitsch /* Start left look right. Looking for e.g. B[0] in A or mergelo. l inclusive, r inclusive. Returns first m such that arr[m] > 126676f2a66SJacob Faibussowitsch x. Output also inclusive. 127676f2a66SJacob Faibussowitsch 128676f2a66SJacob Faibussowitsch NOTE: Both gallopsearch functions CANNOT distinguish between inserting AFTER the entire array vs inserting at entry n!! For example for an array: 129676f2a66SJacob Faibussowitsch 130676f2a66SJacob Faibussowitsch {0,1,2,3,4,5} 131676f2a66SJacob Faibussowitsch 132676f2a66SJacob Faibussowitsch when looking to insert "5" this routine will return *m = 6, but when looking to insert "6" it will ALSO return *m = 6. 133676f2a66SJacob Faibussowitsch */ 134d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscGallopSearchLeft_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m) 135d71ae5a4SJacob Faibussowitsch { 136676f2a66SJacob Faibussowitsch PetscInt last = l, k = 1, mid, cur = l + 1; 137676f2a66SJacob Faibussowitsch 138676f2a66SJacob Faibussowitsch PetscFunctionBegin; 139676f2a66SJacob Faibussowitsch *m = l; 1406bdcaf15SBarry Smith PetscAssert(r >= l, PETSC_COMM_SELF, PETSC_ERR_PLIB, "r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchLeft", r, l); 1419371c9d4SSatish Balay if ((*cmp)(x, arr + r * size, ctx) >= 0) { 1429371c9d4SSatish Balay *m = r; 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1449371c9d4SSatish Balay } 1453ba16761SJacob Faibussowitsch if ((*cmp)(x, (arr) + l * size, ctx) < 0 || PetscUnlikely(!(r - l))) PetscFunctionReturn(PETSC_SUCCESS); 146676f2a66SJacob Faibussowitsch while (PETSC_TRUE) { 1479371c9d4SSatish Balay if (PetscUnlikely(cur > r)) { 1489371c9d4SSatish Balay cur = r; 1499371c9d4SSatish Balay break; 1509371c9d4SSatish Balay } 1514d3610e3SJacob Faibussowitsch if ((*cmp)(x, (arr) + cur * size, ctx) < 0) break; 152676f2a66SJacob Faibussowitsch last = cur; 153676f2a66SJacob Faibussowitsch cur += (k <<= 1) + 1; 154676f2a66SJacob Faibussowitsch ++k; 155676f2a66SJacob Faibussowitsch } 156676f2a66SJacob Faibussowitsch /* standard binary search but take last 0 mid 0 cur 1 into account*/ 157676f2a66SJacob Faibussowitsch while (cur > last + 1) { 158676f2a66SJacob Faibussowitsch mid = last + ((cur - last) >> 1); 1594d3610e3SJacob Faibussowitsch if ((*cmp)(x, (arr) + mid * size, ctx) < 0) { 160676f2a66SJacob Faibussowitsch cur = mid; 161676f2a66SJacob Faibussowitsch } else { 162676f2a66SJacob Faibussowitsch last = mid; 163676f2a66SJacob Faibussowitsch } 164676f2a66SJacob Faibussowitsch } 165676f2a66SJacob Faibussowitsch *m = cur; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167676f2a66SJacob Faibussowitsch } 168676f2a66SJacob Faibussowitsch 169676f2a66SJacob Faibussowitsch /* Start right look left. Looking for e.g. A[-1] in B or mergehi. l inclusive, r inclusive. Returns last m such that arr[m] 170676f2a66SJacob Faibussowitsch < x. Output also inclusive */ 171d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscGallopSearchRight_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m) 172d71ae5a4SJacob Faibussowitsch { 173676f2a66SJacob Faibussowitsch PetscInt last = r, k = 1, mid, cur = r - 1; 174676f2a66SJacob Faibussowitsch 175676f2a66SJacob Faibussowitsch PetscFunctionBegin; 176676f2a66SJacob Faibussowitsch *m = r; 1776bdcaf15SBarry Smith PetscAssert(r >= l, PETSC_COMM_SELF, PETSC_ERR_PLIB, "r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchRight", r, l); 1789371c9d4SSatish Balay if ((*cmp)(x, arr + l * size, ctx) <= 0) { 1799371c9d4SSatish Balay *m = l; 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1819371c9d4SSatish Balay } 1823ba16761SJacob Faibussowitsch if ((*cmp)(x, (arr) + r * size, ctx) > 0 || PetscUnlikely(!(r - l))) PetscFunctionReturn(PETSC_SUCCESS); 183676f2a66SJacob Faibussowitsch while (PETSC_TRUE) { 1849371c9d4SSatish Balay if (PetscUnlikely(cur < l)) { 1859371c9d4SSatish Balay cur = l; 1869371c9d4SSatish Balay break; 1879371c9d4SSatish Balay } 1884d3610e3SJacob Faibussowitsch if ((*cmp)(x, (arr) + cur * size, ctx) > 0) break; 189676f2a66SJacob Faibussowitsch last = cur; 190676f2a66SJacob Faibussowitsch cur -= (k <<= 1) + 1; 191676f2a66SJacob Faibussowitsch ++k; 192676f2a66SJacob Faibussowitsch } 193676f2a66SJacob Faibussowitsch /* standard binary search but take last r-1 mid r-1 cur r-2 into account*/ 194676f2a66SJacob Faibussowitsch while (last > cur + 1) { 195676f2a66SJacob Faibussowitsch mid = last - ((last - cur) >> 1); 1964d3610e3SJacob Faibussowitsch if ((*cmp)(x, (arr) + mid * size, ctx) > 0) { 197676f2a66SJacob Faibussowitsch cur = mid; 198676f2a66SJacob Faibussowitsch } else { 199676f2a66SJacob Faibussowitsch last = mid; 200676f2a66SJacob Faibussowitsch } 201676f2a66SJacob Faibussowitsch } 202676f2a66SJacob Faibussowitsch *m = cur; 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204676f2a66SJacob Faibussowitsch } 205676f2a66SJacob Faibussowitsch 206676f2a66SJacob Faibussowitsch /* Mergesort where size of left half <= size of right half, so mergesort is done left to right. Arr should be pointer to 207676f2a66SJacob Faibussowitsch complete array, left is first index of left array, mid is first index of right array, right is last index of right 208676f2a66SJacob Faibussowitsch array */ 209d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeLo_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) 210d71ae5a4SJacob Faibussowitsch { 211676f2a66SJacob Faibussowitsch PetscInt i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0; 212676f2a66SJacob Faibussowitsch const PetscInt llen = mid - left; 213676f2a66SJacob Faibussowitsch 214676f2a66SJacob Faibussowitsch PetscFunctionBegin; 215676f2a66SJacob Faibussowitsch Petsc_memcpy(tarr, arr + (left * size), llen * size); 216676f2a66SJacob Faibussowitsch while ((i < llen) && (j <= right)) { 2174d3610e3SJacob Faibussowitsch if ((*cmp)(tarr + (i * size), arr + (j * size), ctx) < 0) { 218676f2a66SJacob Faibussowitsch Petsc_memcpy(arr + (k * size), tarr + (i * size), size); 219676f2a66SJacob Faibussowitsch ++k; 220676f2a66SJacob Faibussowitsch gallopright = 0; 221676f2a66SJacob Faibussowitsch if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) { 222676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 223676f2a66SJacob Faibussowitsch do { 224676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 225676f2a66SJacob Faibussowitsch /* search temp for right[j], can move up to that of temp into arr immediately */ 2269566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen - 1, arr + (j * size), &l1)); 227676f2a66SJacob Faibussowitsch diff1 = l1 - i; 228676f2a66SJacob Faibussowitsch Petsc_memcpy(arr + (k * size), tarr + (i * size), diff1 * size); 229676f2a66SJacob Faibussowitsch k += diff1; 230676f2a66SJacob Faibussowitsch i = l1; 231676f2a66SJacob Faibussowitsch /* search right for temp[i], can move up to that many of right into arr */ 2329566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr + (i * size), &l2)); 233676f2a66SJacob Faibussowitsch diff2 = l2 - j; 234676f2a66SJacob Faibussowitsch Petsc_memmove((arr) + k * size, (arr) + j * size, diff2 * size); 235676f2a66SJacob Faibussowitsch k += diff2; 236676f2a66SJacob Faibussowitsch j = l2; 237676f2a66SJacob Faibussowitsch if (i >= llen || j > right) break; 238676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 239676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 240676f2a66SJacob Faibussowitsch } 241676f2a66SJacob Faibussowitsch } else { 242676f2a66SJacob Faibussowitsch Petsc_memmove(arr + (k * size), arr + (j * size), size); 243676f2a66SJacob Faibussowitsch ++k; 244676f2a66SJacob Faibussowitsch gallopleft = 0; 245676f2a66SJacob Faibussowitsch if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) { 246676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 247676f2a66SJacob Faibussowitsch do { 248676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 249676f2a66SJacob Faibussowitsch /* search right for temp[i], can move up to that many of right into arr */ 2509566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr + (i * size), &l2)); 251676f2a66SJacob Faibussowitsch diff2 = l2 - j; 252676f2a66SJacob Faibussowitsch Petsc_memmove(arr + (k * size), arr + (j * size), diff2 * size); 253676f2a66SJacob Faibussowitsch k += diff2; 254676f2a66SJacob Faibussowitsch j = l2; 255676f2a66SJacob Faibussowitsch /* search temp for right[j], can copy up to that of temp into arr immediately */ 2569566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen - 1, arr + (j * size), &l1)); 257676f2a66SJacob Faibussowitsch diff1 = l1 - i; 258676f2a66SJacob Faibussowitsch Petsc_memcpy(arr + (k * size), tarr + (i * size), diff1 * size); 259676f2a66SJacob Faibussowitsch k += diff1; 260676f2a66SJacob Faibussowitsch i = l1; 261676f2a66SJacob Faibussowitsch if (i >= llen || j > right) break; 262676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 263676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 264676f2a66SJacob Faibussowitsch } 265676f2a66SJacob Faibussowitsch } 266676f2a66SJacob Faibussowitsch } 267ad540459SPierre Jolivet if (i < llen) Petsc_memcpy(arr + (k * size), tarr + (i * size), (llen - i) * size); 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269676f2a66SJacob Faibussowitsch } 270676f2a66SJacob Faibussowitsch 271d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeLoWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) 272d71ae5a4SJacob Faibussowitsch { 273676f2a66SJacob Faibussowitsch PetscInt i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0; 274676f2a66SJacob Faibussowitsch const PetscInt llen = mid - left; 275676f2a66SJacob Faibussowitsch 276676f2a66SJacob Faibussowitsch PetscFunctionBegin; 277676f2a66SJacob Faibussowitsch Petsc_memcpy2(atarr, arr + (left * asize), llen * asize, btarr, barr + (left * bsize), llen * bsize); 278676f2a66SJacob Faibussowitsch while ((i < llen) && (j <= right)) { 2794d3610e3SJacob Faibussowitsch if ((*cmp)(atarr + (i * asize), arr + (j * asize), ctx) < 0) { 280676f2a66SJacob Faibussowitsch Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), asize, barr + (k * bsize), btarr + (i * bsize), bsize); 281676f2a66SJacob Faibussowitsch ++k; 282676f2a66SJacob Faibussowitsch gallopright = 0; 283676f2a66SJacob Faibussowitsch if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) { 284676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 285676f2a66SJacob Faibussowitsch do { 286676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 287676f2a66SJacob Faibussowitsch /* search temp for right[j], can move up to that of temp into arr immediately */ 2889566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen - 1, arr + (j * asize), &l1)); 289676f2a66SJacob Faibussowitsch diff1 = l1 - i; 290676f2a66SJacob Faibussowitsch Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), diff1 * asize, barr + (k * bsize), btarr + (i * bsize), diff1 * bsize); 291676f2a66SJacob Faibussowitsch k += diff1; 292676f2a66SJacob Faibussowitsch i = l1; 293676f2a66SJacob Faibussowitsch /* search right for temp[i], can move up to that many of right into arr */ 2949566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr + (i * asize), &l2)); 295676f2a66SJacob Faibussowitsch diff2 = l2 - j; 296676f2a66SJacob Faibussowitsch Petsc_memmove2(arr + (k * asize), arr + (j * asize), diff2 * asize, barr + (k * bsize), barr + (j * bsize), diff2 * bsize); 297676f2a66SJacob Faibussowitsch k += diff2; 298676f2a66SJacob Faibussowitsch j = l2; 299676f2a66SJacob Faibussowitsch if (i >= llen || j > right) break; 300676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 301676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 302676f2a66SJacob Faibussowitsch } 303676f2a66SJacob Faibussowitsch } else { 304676f2a66SJacob Faibussowitsch Petsc_memmove2(arr + (k * asize), arr + (j * asize), asize, barr + (k * bsize), barr + (j * bsize), bsize); 305676f2a66SJacob Faibussowitsch ++k; 306676f2a66SJacob Faibussowitsch gallopleft = 0; 307676f2a66SJacob Faibussowitsch if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) { 308676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 309676f2a66SJacob Faibussowitsch do { 310676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 311676f2a66SJacob Faibussowitsch /* search right for temp[i], can move up to that many of right into arr */ 3129566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr + (i * asize), &l2)); 313676f2a66SJacob Faibussowitsch diff2 = l2 - j; 314676f2a66SJacob Faibussowitsch Petsc_memmove2(arr + (k * asize), arr + (j * asize), diff2 * asize, barr + (k * bsize), barr + (j * bsize), diff2 * bsize); 315676f2a66SJacob Faibussowitsch k += diff2; 316676f2a66SJacob Faibussowitsch j = l2; 317676f2a66SJacob Faibussowitsch /* search temp for right[j], can copy up to that of temp into arr immediately */ 3189566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen - 1, arr + (j * asize), &l1)); 319676f2a66SJacob Faibussowitsch diff1 = l1 - i; 320676f2a66SJacob Faibussowitsch Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), diff1 * asize, barr + (k * bsize), btarr + (i * bsize), diff1 * bsize); 321676f2a66SJacob Faibussowitsch k += diff1; 322676f2a66SJacob Faibussowitsch i = l1; 323676f2a66SJacob Faibussowitsch if (i >= llen || j > right) break; 324676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 325676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 326676f2a66SJacob Faibussowitsch } 327676f2a66SJacob Faibussowitsch } 328676f2a66SJacob Faibussowitsch } 329ad540459SPierre Jolivet if (i < llen) Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), (llen - i) * asize, barr + (k * bsize), btarr + (i * bsize), (llen - i) * bsize); 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 331676f2a66SJacob Faibussowitsch } 332676f2a66SJacob Faibussowitsch 333676f2a66SJacob Faibussowitsch /* Mergesort where size of right half < size of left half, so mergesort is done right to left. Arr should be pointer to 334676f2a66SJacob Faibussowitsch complete array, left is first index of left array, mid is first index of right array, right is last index of right 335676f2a66SJacob Faibussowitsch array */ 336d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeHi_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) 337d71ae5a4SJacob Faibussowitsch { 338676f2a66SJacob Faibussowitsch PetscInt i = right - mid, j = mid - 1, k = right, gallopleft = 0, gallopright = 0; 339676f2a66SJacob Faibussowitsch const PetscInt rlen = right - mid + 1; 340676f2a66SJacob Faibussowitsch 341676f2a66SJacob Faibussowitsch PetscFunctionBegin; 342676f2a66SJacob Faibussowitsch Petsc_memcpy(tarr, (arr) + mid * size, rlen * size); 343676f2a66SJacob Faibussowitsch while ((i >= 0) && (j >= left)) { 3444d3610e3SJacob Faibussowitsch if ((*cmp)((tarr) + i * size, (arr) + j * size, ctx) > 0) { 345676f2a66SJacob Faibussowitsch Petsc_memcpy((arr) + k * size, (tarr) + i * size, size); 346676f2a66SJacob Faibussowitsch --k; 347676f2a66SJacob Faibussowitsch gallopleft = 0; 348676f2a66SJacob Faibussowitsch if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) { 349676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 350676f2a66SJacob Faibussowitsch do { 351676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 352676f2a66SJacob Faibussowitsch /* search temp for left[j], can copy up to that many of temp into arr */ 3539566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr) + j * size, &l1)); 354676f2a66SJacob Faibussowitsch diff1 = i - l1; 355676f2a66SJacob Faibussowitsch Petsc_memcpy((arr) + (k - diff1 + 1) * size, (tarr) + (l1 + 1) * size, diff1 * size); 356676f2a66SJacob Faibussowitsch k -= diff1; 357676f2a66SJacob Faibussowitsch i = l1; 358676f2a66SJacob Faibussowitsch /* search left for temp[i], can move up to that many of left up arr */ 3599566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr) + i * size, &l2)); 360676f2a66SJacob Faibussowitsch diff2 = j - l2; 361676f2a66SJacob Faibussowitsch Petsc_memmove((arr) + (k - diff2 + 1) * size, (arr) + (l2 + 1) * size, diff2 * size); 362676f2a66SJacob Faibussowitsch k -= diff2; 363676f2a66SJacob Faibussowitsch j = l2; 364676f2a66SJacob Faibussowitsch if (i < 0 || j < left) break; 365676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 366676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 367676f2a66SJacob Faibussowitsch } 368676f2a66SJacob Faibussowitsch } else { 369676f2a66SJacob Faibussowitsch Petsc_memmove((arr) + k * size, (arr) + j * size, size); 370676f2a66SJacob Faibussowitsch --k; 371676f2a66SJacob Faibussowitsch gallopright = 0; 372676f2a66SJacob Faibussowitsch if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) { 373676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 374676f2a66SJacob Faibussowitsch do { 375676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 376676f2a66SJacob Faibussowitsch /* search left for temp[i], can move up to that many of left up arr */ 3779566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr) + i * size, &l2)); 378676f2a66SJacob Faibussowitsch diff2 = j - l2; 379676f2a66SJacob Faibussowitsch Petsc_memmove((arr) + (k - diff2 + 1) * size, (arr) + (l2 + 1) * size, diff2 * size); 380676f2a66SJacob Faibussowitsch k -= diff2; 381676f2a66SJacob Faibussowitsch j = l2; 382676f2a66SJacob Faibussowitsch /* search temp for left[j], can copy up to that many of temp into arr */ 3839566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr) + j * size, &l1)); 384676f2a66SJacob Faibussowitsch diff1 = i - l1; 385676f2a66SJacob Faibussowitsch Petsc_memcpy((arr) + (k - diff1 + 1) * size, (tarr) + (l1 + 1) * size, diff1 * size); 386676f2a66SJacob Faibussowitsch k -= diff1; 387676f2a66SJacob Faibussowitsch i = l1; 388676f2a66SJacob Faibussowitsch if (i < 0 || j < left) break; 389676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 390676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 391676f2a66SJacob Faibussowitsch } 392676f2a66SJacob Faibussowitsch } 393676f2a66SJacob Faibussowitsch } 394ad540459SPierre Jolivet if (i >= 0) Petsc_memcpy((arr) + left * size, tarr, (i + 1) * size); 3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 396676f2a66SJacob Faibussowitsch } 397676f2a66SJacob Faibussowitsch 398d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeHiWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) 399d71ae5a4SJacob Faibussowitsch { 400676f2a66SJacob Faibussowitsch PetscInt i = right - mid, j = mid - 1, k = right, gallopleft = 0, gallopright = 0; 401676f2a66SJacob Faibussowitsch const PetscInt rlen = right - mid + 1; 402676f2a66SJacob Faibussowitsch 403676f2a66SJacob Faibussowitsch PetscFunctionBegin; 404676f2a66SJacob Faibussowitsch Petsc_memcpy2(atarr, (arr) + mid * asize, rlen * asize, btarr, (barr) + mid * bsize, rlen * bsize); 405676f2a66SJacob Faibussowitsch while ((i >= 0) && (j >= left)) { 4064d3610e3SJacob Faibussowitsch if ((*cmp)((atarr) + i * asize, (arr) + j * asize, ctx) > 0) { 407676f2a66SJacob Faibussowitsch Petsc_memcpy2((arr) + k * asize, (atarr) + i * asize, asize, (barr) + k * bsize, (btarr) + i * bsize, bsize); 408676f2a66SJacob Faibussowitsch --k; 409676f2a66SJacob Faibussowitsch gallopleft = 0; 410676f2a66SJacob Faibussowitsch if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) { 411676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 412676f2a66SJacob Faibussowitsch do { 413676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 414676f2a66SJacob Faibussowitsch /* search temp for left[j], can copy up to that many of temp into arr */ 4159566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr) + j * asize, &l1)); 416676f2a66SJacob Faibussowitsch diff1 = i - l1; 417676f2a66SJacob Faibussowitsch Petsc_memcpy2((arr) + (k - diff1 + 1) * asize, (atarr) + (l1 + 1) * asize, diff1 * asize, (barr) + (k - diff1 + 1) * bsize, (btarr) + (l1 + 1) * bsize, diff1 * bsize); 418676f2a66SJacob Faibussowitsch k -= diff1; 419676f2a66SJacob Faibussowitsch i = l1; 420676f2a66SJacob Faibussowitsch /* search left for temp[i], can move up to that many of left up arr */ 4219566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr) + i * asize, &l2)); 422676f2a66SJacob Faibussowitsch diff2 = j - l2; 423676f2a66SJacob Faibussowitsch Petsc_memmove2((arr) + (k - diff2 + 1) * asize, (arr) + (l2 + 1) * asize, diff2 * asize, (barr) + (k - diff2 + 1) * bsize, (barr) + (l2 + 1) * bsize, diff2 * bsize); 424676f2a66SJacob Faibussowitsch k -= diff2; 425676f2a66SJacob Faibussowitsch j = l2; 426676f2a66SJacob Faibussowitsch if (i < 0 || j < left) break; 427676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 428676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 429676f2a66SJacob Faibussowitsch } 430676f2a66SJacob Faibussowitsch } else { 431676f2a66SJacob Faibussowitsch Petsc_memmove2((arr) + k * asize, (arr) + j * asize, asize, (barr) + k * bsize, (barr) + j * bsize, bsize); 432676f2a66SJacob Faibussowitsch --k; 433676f2a66SJacob Faibussowitsch gallopright = 0; 434676f2a66SJacob Faibussowitsch if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) { 435676f2a66SJacob Faibussowitsch PetscInt l1, l2, diff1, diff2; 436676f2a66SJacob Faibussowitsch do { 437676f2a66SJacob Faibussowitsch if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL; 438676f2a66SJacob Faibussowitsch /* search left for temp[i], can move up to that many of left up arr */ 4399566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr) + i * asize, &l2)); 440676f2a66SJacob Faibussowitsch diff2 = j - l2; 441676f2a66SJacob Faibussowitsch Petsc_memmove2((arr) + (k - diff2 + 1) * asize, (arr) + (l2 + 1) * asize, diff2 * asize, (barr) + (k - diff2 + 1) * bsize, (barr) + (l2 + 1) * bsize, diff2 * bsize); 442676f2a66SJacob Faibussowitsch k -= diff2; 443676f2a66SJacob Faibussowitsch j = l2; 444676f2a66SJacob Faibussowitsch /* search temp for left[j], can copy up to that many of temp into arr */ 4459566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr) + j * asize, &l1)); 446676f2a66SJacob Faibussowitsch diff1 = i - l1; 447676f2a66SJacob Faibussowitsch Petsc_memcpy2((arr) + (k - diff1 + 1) * asize, (atarr) + (l1 + 1) * asize, diff1 * asize, (barr) + (k - diff1 + 1) * bsize, (btarr) + (l1 + 1) * bsize, diff1 * bsize); 448676f2a66SJacob Faibussowitsch k -= diff1; 449676f2a66SJacob Faibussowitsch i = l1; 450676f2a66SJacob Faibussowitsch if (i < 0 || j < left) break; 451676f2a66SJacob Faibussowitsch } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL); 452676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; 453676f2a66SJacob Faibussowitsch } 454676f2a66SJacob Faibussowitsch } 455676f2a66SJacob Faibussowitsch } 456ad540459SPierre Jolivet if (i >= 0) Petsc_memcpy2((arr) + left * asize, atarr, (i + 1) * asize, (barr) + left * bsize, btarr, (i + 1) * bsize); 4573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 458676f2a66SJacob Faibussowitsch } 459676f2a66SJacob Faibussowitsch 460676f2a66SJacob Faibussowitsch /* Left is inclusive lower bound of array slice, start is start location of unsorted section, right is inclusive upper 461676f2a66SJacob Faibussowitsch bound of array slice. If unsure of where unsorted section starts or if entire length is unsorted pass start = left */ 462d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) 463d71ae5a4SJacob Faibussowitsch { 464676f2a66SJacob Faibussowitsch PetscInt i = start == left ? start + 1 : start; 465676f2a66SJacob Faibussowitsch 466676f2a66SJacob Faibussowitsch PetscFunctionBegin; 467676f2a66SJacob Faibussowitsch for (; i <= right; ++i) { 468676f2a66SJacob Faibussowitsch PetscInt j = i - 1; 469676f2a66SJacob Faibussowitsch Petsc_memcpy(tarr, arr + (i * size), size); 4704d3610e3SJacob Faibussowitsch while ((j >= left) && ((*cmp)(tarr, (arr) + j * size, ctx) < 0)) { 471676f2a66SJacob Faibussowitsch COPYSWAPPY(arr + (j + 1) * size, arr + j * size, tarr + size, size); 472676f2a66SJacob Faibussowitsch --j; 473676f2a66SJacob Faibussowitsch } 474676f2a66SJacob Faibussowitsch Petsc_memcpy((arr) + (j + 1) * size, tarr, size); 475676f2a66SJacob Faibussowitsch } 4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 477676f2a66SJacob Faibussowitsch } 478676f2a66SJacob Faibussowitsch 479d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) 480d71ae5a4SJacob Faibussowitsch { 481676f2a66SJacob Faibussowitsch PetscInt i = start == left ? start + 1 : start; 482676f2a66SJacob Faibussowitsch 483676f2a66SJacob Faibussowitsch PetscFunctionBegin; 484676f2a66SJacob Faibussowitsch for (; i <= right; ++i) { 485676f2a66SJacob Faibussowitsch PetscInt j = i - 1; 486676f2a66SJacob Faibussowitsch Petsc_memcpy2(atarr, arr + (i * asize), asize, btarr, barr + (i * bsize), bsize); 4874d3610e3SJacob Faibussowitsch while ((j >= left) && ((*cmp)(atarr, arr + (j * asize), ctx) < 0)) { 488676f2a66SJacob Faibussowitsch COPYSWAPPY2(arr + (j + 1) * asize, arr + j * asize, asize, barr + (j + 1) * bsize, barr + j * bsize, bsize, atarr + asize); 489676f2a66SJacob Faibussowitsch --j; 490676f2a66SJacob Faibussowitsch } 491676f2a66SJacob Faibussowitsch Petsc_memcpy2(arr + (j + 1) * asize, atarr, asize, barr + (j + 1) * bsize, btarr, bsize); 492676f2a66SJacob Faibussowitsch } 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 494676f2a66SJacob Faibussowitsch } 495676f2a66SJacob Faibussowitsch 496676f2a66SJacob Faibussowitsch /* See PetscInsertionSort_Private */ 497d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBinaryInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) 498d71ae5a4SJacob Faibussowitsch { 499676f2a66SJacob Faibussowitsch PetscInt i = start == left ? start + 1 : start; 500676f2a66SJacob Faibussowitsch 501676f2a66SJacob Faibussowitsch PetscFunctionBegin; 502676f2a66SJacob Faibussowitsch for (; i <= right; ++i) { 503676f2a66SJacob Faibussowitsch PetscInt l = left, r = i; 504676f2a66SJacob Faibussowitsch Petsc_memcpy(tarr, arr + (i * size), size); 505676f2a66SJacob Faibussowitsch do { 506676f2a66SJacob Faibussowitsch const PetscInt m = l + ((r - l) >> 1); 5074d3610e3SJacob Faibussowitsch if ((*cmp)(tarr, arr + (m * size), ctx) < 0) { 508676f2a66SJacob Faibussowitsch r = m; 509676f2a66SJacob Faibussowitsch } else { 510676f2a66SJacob Faibussowitsch l = m + 1; 511676f2a66SJacob Faibussowitsch } 512676f2a66SJacob Faibussowitsch } while (l < r); 513676f2a66SJacob Faibussowitsch Petsc_memmove(arr + ((l + 1) * size), arr + (l * size), (i - l) * size); 514676f2a66SJacob Faibussowitsch Petsc_memcpy(arr + (l * size), tarr, size); 515676f2a66SJacob Faibussowitsch } 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 517676f2a66SJacob Faibussowitsch } 518676f2a66SJacob Faibussowitsch 519d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBinaryInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) 520d71ae5a4SJacob Faibussowitsch { 521676f2a66SJacob Faibussowitsch PetscInt i = start == left ? start + 1 : start; 522676f2a66SJacob Faibussowitsch 523676f2a66SJacob Faibussowitsch PetscFunctionBegin; 524676f2a66SJacob Faibussowitsch for (; i <= right; ++i) { 525676f2a66SJacob Faibussowitsch PetscInt l = left, r = i; 526676f2a66SJacob Faibussowitsch Petsc_memcpy2(atarr, arr + (i * asize), asize, btarr, barr + (i * bsize), bsize); 527676f2a66SJacob Faibussowitsch do { 528676f2a66SJacob Faibussowitsch const PetscInt m = l + ((r - l) >> 1); 5294d3610e3SJacob Faibussowitsch if ((*cmp)(atarr, arr + (m * asize), ctx) < 0) { 530676f2a66SJacob Faibussowitsch r = m; 531676f2a66SJacob Faibussowitsch } else { 532676f2a66SJacob Faibussowitsch l = m + 1; 533676f2a66SJacob Faibussowitsch } 534676f2a66SJacob Faibussowitsch } while (l < r); 535676f2a66SJacob Faibussowitsch Petsc_memmove2(arr + ((l + 1) * asize), arr + (l * asize), (i - l) * asize, barr + ((l + 1) * bsize), barr + (l * bsize), (i - l) * bsize); 536676f2a66SJacob Faibussowitsch Petsc_memcpy2(arr + (l * asize), atarr, asize, barr + (l * bsize), btarr, bsize); 537676f2a66SJacob Faibussowitsch } 5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 539676f2a66SJacob Faibussowitsch } 540676f2a66SJacob Faibussowitsch 541676f2a66SJacob Faibussowitsch typedef struct { 542676f2a66SJacob Faibussowitsch PetscInt size; 543676f2a66SJacob Faibussowitsch PetscInt start; 54489d292beSStefano Zampini #if defined(__CRAY_AARCH64) /* segfaults with Cray compilers for aarch64 on a64FX */ 54589d292beSStefano Zampini } PetscTimSortStack; 54689d292beSStefano Zampini #else 547676f2a66SJacob Faibussowitsch } PetscTimSortStack PETSC_ATTRIBUTEALIGNED(2 * sizeof(PetscInt)); 54889d292beSStefano Zampini #endif 549676f2a66SJacob Faibussowitsch 550676f2a66SJacob Faibussowitsch typedef struct { 551676f2a66SJacob Faibussowitsch char *ptr PETSC_ATTRIBUTEALIGNED(PETSC_MEMALIGN); 552676f2a66SJacob Faibussowitsch size_t size; 553676f2a66SJacob Faibussowitsch size_t maxsize; 554676f2a66SJacob Faibussowitsch } PetscTimSortBuffer; 555676f2a66SJacob Faibussowitsch 556d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortResizeBuffer_Private(PetscTimSortBuffer *buff, size_t newSize) 557d71ae5a4SJacob Faibussowitsch { 558676f2a66SJacob Faibussowitsch PetscFunctionBegin; 5593ba16761SJacob Faibussowitsch if (PetscLikely(newSize <= buff->size)) PetscFunctionReturn(PETSC_SUCCESS); 560676f2a66SJacob Faibussowitsch { 561676f2a66SJacob Faibussowitsch /* Can't be larger than n, there is merit to simply allocating buff to n to begin with */ 562676f2a66SJacob Faibussowitsch size_t newMax = PetscMin(newSize * newSize, buff->maxsize); 5639566063dSJacob Faibussowitsch PetscCall(PetscFree(buff->ptr)); 5649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(newMax, &buff->ptr)); 565676f2a66SJacob Faibussowitsch buff->size = newMax; 566676f2a66SJacob Faibussowitsch } 5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 568676f2a66SJacob Faibussowitsch } 569676f2a66SJacob Faibussowitsch 570d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortForceCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt stacksize) 571d71ae5a4SJacob Faibussowitsch { 572676f2a66SJacob Faibussowitsch PetscFunctionBegin; 573676f2a66SJacob Faibussowitsch for (; stacksize; --stacksize) { 574676f2a66SJacob Faibussowitsch /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */ 5754d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[stacksize].start - 1) * size, arr + (stack[stacksize].start) * size, ctx) > 0) { 576676f2a66SJacob Faibussowitsch PetscInt l, m = stack[stacksize].start, r; 577676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 5789566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[stacksize - 1].start, stack[stacksize].start - 1, (arr) + (stack[stacksize].start) * size, &l)); 579676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 5809566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[stacksize].start, stack[stacksize].start + stack[stacksize].size - 1, (arr) + (stack[stacksize].start - 1) * size, &r)); 581676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 5829566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size)); 5839566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 584676f2a66SJacob Faibussowitsch } else { 5859566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size)); 5869566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 587676f2a66SJacob Faibussowitsch } 588676f2a66SJacob Faibussowitsch } 589676f2a66SJacob Faibussowitsch /* Update A with merge */ 590676f2a66SJacob Faibussowitsch stack[stacksize - 1].size += stack[stacksize].size; 591676f2a66SJacob Faibussowitsch } 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593676f2a66SJacob Faibussowitsch } 594676f2a66SJacob Faibussowitsch 595d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortForceCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, void *ctx, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt stacksize) 596d71ae5a4SJacob Faibussowitsch { 597676f2a66SJacob Faibussowitsch PetscFunctionBegin; 598676f2a66SJacob Faibussowitsch for (; stacksize; --stacksize) { 599676f2a66SJacob Faibussowitsch /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */ 6004d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[stacksize].start - 1) * asize, arr + (stack[stacksize].start) * asize, ctx) > 0) { 601676f2a66SJacob Faibussowitsch PetscInt l, m = stack[stacksize].start, r; 602676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 6039566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[stacksize - 1].start, stack[stacksize].start - 1, (arr) + (stack[stacksize].start) * asize, &l)); 604676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 6059566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[stacksize].start, stack[stacksize].start + stack[stacksize].size - 1, (arr) + (stack[stacksize].start - 1) * asize, &r)); 606676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 6079566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize)); 6089566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize)); 6099566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 610676f2a66SJacob Faibussowitsch } else { 6119566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize)); 6129566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize)); 6139566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 614676f2a66SJacob Faibussowitsch } 615676f2a66SJacob Faibussowitsch } 616676f2a66SJacob Faibussowitsch /* Update A with merge */ 617676f2a66SJacob Faibussowitsch stack[stacksize - 1].size += stack[stacksize].size; 618676f2a66SJacob Faibussowitsch } 6193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 620676f2a66SJacob Faibussowitsch } 621676f2a66SJacob Faibussowitsch 622d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt *stacksize) 623d71ae5a4SJacob Faibussowitsch { 624676f2a66SJacob Faibussowitsch PetscInt i = *stacksize; 625676f2a66SJacob Faibussowitsch 626676f2a66SJacob Faibussowitsch PetscFunctionBegin; 627676f2a66SJacob Faibussowitsch while (i) { 628676f2a66SJacob Faibussowitsch PetscInt l, m, r, itemp = i; 629676f2a66SJacob Faibussowitsch 630676f2a66SJacob Faibussowitsch if (i == 1) { 631676f2a66SJacob Faibussowitsch /* A = stack[i-1], B = stack[i] */ 632676f2a66SJacob Faibussowitsch if (stack[i - 1].size < stack[i].size) { 633676f2a66SJacob Faibussowitsch /* if A[-1] <= B[0] then sorted */ 6344d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i].start - 1) * size, arr + (stack[i].start) * size, ctx) > 0) { 635676f2a66SJacob Faibussowitsch m = stack[i].start; 636676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 6379566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 1].start, stack[i].start - 1, (arr) + stack[i].start * size, &l)); 638676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 6399566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[i].start, stack[i].start + stack[i].size - 1, arr + (stack[i].start - 1) * size, &r)); 640676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 6419566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size)); 6429566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 643676f2a66SJacob Faibussowitsch } else { 6449566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size)); 6459566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 646676f2a66SJacob Faibussowitsch } 647676f2a66SJacob Faibussowitsch } 648676f2a66SJacob Faibussowitsch /* Update A with merge */ 649676f2a66SJacob Faibussowitsch stack[i - 1].size += stack[i].size; 650676f2a66SJacob Faibussowitsch --i; 651676f2a66SJacob Faibussowitsch } 652676f2a66SJacob Faibussowitsch } else { 653676f2a66SJacob Faibussowitsch /* i > 2, i.e. C exists 654676f2a66SJacob Faibussowitsch A = stack[i-2], B = stack[i-1], C = stack[i]; */ 655676f2a66SJacob Faibussowitsch if (stack[i - 2].size <= stack[i - 1].size + stack[i].size) { 656676f2a66SJacob Faibussowitsch if (stack[i - 2].size < stack[i].size) { 657676f2a66SJacob Faibussowitsch /* merge B into A, but if A[-1] <= B[0] then already sorted */ 6584d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i - 1].start - 1) * size, arr + (stack[i - 1].start) * size, ctx) > 0) { 659676f2a66SJacob Faibussowitsch m = stack[i - 1].start; 660676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 6619566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 2].start, stack[i - 1].start - 1, (arr) + (stack[i - 1].start) * size, &l)); 662676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 6639566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[i - 1].start, stack[i - 1].start + stack[i - 1].size - 1, (arr) + (stack[i - 1].start - 1) * size, &r)); 664676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 6659566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size)); 6669566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 667676f2a66SJacob Faibussowitsch } else { 6689566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size)); 6699566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 670676f2a66SJacob Faibussowitsch } 671676f2a66SJacob Faibussowitsch } 672676f2a66SJacob Faibussowitsch /* Update A with merge */ 673676f2a66SJacob Faibussowitsch stack[i - 2].size += stack[i - 1].size; 674676f2a66SJacob Faibussowitsch /* Push C up the stack */ 675676f2a66SJacob Faibussowitsch stack[i - 1].start = stack[i].start; 676676f2a66SJacob Faibussowitsch stack[i - 1].size = stack[i].size; 677676f2a66SJacob Faibussowitsch } else { 678676f2a66SJacob Faibussowitsch /* merge C into B */ 679676f2a66SJacob Faibussowitsch mergeBC: 680676f2a66SJacob Faibussowitsch /* If B[-1] <= C[0] then... you know the drill */ 6814d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i].start - 1) * size, arr + (stack[i].start) * size, ctx) > 0) { 682676f2a66SJacob Faibussowitsch m = stack[i].start; 683676f2a66SJacob Faibussowitsch /* Search B for C[0] insertion */ 6849566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 1].start, stack[i].start - 1, arr + stack[i].start * size, &l)); 685676f2a66SJacob Faibussowitsch /* Search C for B[-1] insertion */ 6869566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[i].start, stack[i].start + stack[i].size - 1, (arr) + (stack[i].start - 1) * size, &r)); 687676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 6889566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size)); 6899566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 690676f2a66SJacob Faibussowitsch } else { 6919566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size)); 6929566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r)); 693676f2a66SJacob Faibussowitsch } 694676f2a66SJacob Faibussowitsch } 695676f2a66SJacob Faibussowitsch /* Update B with merge */ 696676f2a66SJacob Faibussowitsch stack[i - 1].size += stack[i].size; 697676f2a66SJacob Faibussowitsch } 698676f2a66SJacob Faibussowitsch --i; 699676f2a66SJacob Faibussowitsch } else if (stack[i - 1].size <= stack[i].size) { 700676f2a66SJacob Faibussowitsch /* merge C into B */ 701676f2a66SJacob Faibussowitsch goto mergeBC; 702676f2a66SJacob Faibussowitsch } 703676f2a66SJacob Faibussowitsch } 704676f2a66SJacob Faibussowitsch if (itemp == i) break; 705676f2a66SJacob Faibussowitsch } 706676f2a66SJacob Faibussowitsch *stacksize = i; 7073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 708676f2a66SJacob Faibussowitsch } 709676f2a66SJacob Faibussowitsch 710d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, void *ctx, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt *stacksize) 711d71ae5a4SJacob Faibussowitsch { 712676f2a66SJacob Faibussowitsch PetscInt i = *stacksize; 713676f2a66SJacob Faibussowitsch 714676f2a66SJacob Faibussowitsch PetscFunctionBegin; 715676f2a66SJacob Faibussowitsch while (i) { 716676f2a66SJacob Faibussowitsch PetscInt l, m, r, itemp = i; 717676f2a66SJacob Faibussowitsch 718676f2a66SJacob Faibussowitsch if (i == 1) { 719676f2a66SJacob Faibussowitsch /* A = stack[i-1], B = stack[i] */ 720676f2a66SJacob Faibussowitsch if (stack[i - 1].size < stack[i].size) { 721676f2a66SJacob Faibussowitsch /* if A[-1] <= B[0] then sorted */ 7224d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i].start - 1) * asize, arr + (stack[i].start) * asize, ctx) > 0) { 723676f2a66SJacob Faibussowitsch m = stack[i].start; 724676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 7259566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 1].start, stack[i].start - 1, (arr) + stack[i].start * asize, &l)); 726676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 7279566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[i].start, stack[i].start + stack[i].size - 1, arr + (stack[i].start - 1) * asize, &r)); 728676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 7299566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize)); 7309566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize)); 7319566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 732676f2a66SJacob Faibussowitsch } else { 7339566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize)); 7349566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize)); 7359566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 736676f2a66SJacob Faibussowitsch } 737676f2a66SJacob Faibussowitsch } 738676f2a66SJacob Faibussowitsch /* Update A with merge */ 739676f2a66SJacob Faibussowitsch stack[i - 1].size += stack[i].size; 740676f2a66SJacob Faibussowitsch --i; 741676f2a66SJacob Faibussowitsch } 742676f2a66SJacob Faibussowitsch } else { 743676f2a66SJacob Faibussowitsch /* i > 2, i.e. C exists 744676f2a66SJacob Faibussowitsch A = stack[i-2], B = stack[i-1], C = stack[i]; */ 745676f2a66SJacob Faibussowitsch if (stack[i - 2].size <= stack[i - 1].size + stack[i].size) { 746676f2a66SJacob Faibussowitsch if (stack[i - 2].size < stack[i].size) { 747676f2a66SJacob Faibussowitsch /* merge B into A, but if A[-1] <= B[0] then already sorted */ 7484d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i - 1].start - 1) * asize, arr + (stack[i - 1].start) * asize, ctx) > 0) { 749676f2a66SJacob Faibussowitsch m = stack[i - 1].start; 750676f2a66SJacob Faibussowitsch /* Search A for B[0] insertion */ 7519566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 2].start, stack[i - 1].start - 1, (arr) + (stack[i - 1].start) * asize, &l)); 752676f2a66SJacob Faibussowitsch /* Search B for A[-1] insertion */ 7539566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[i - 1].start, stack[i - 1].start + stack[i - 1].size - 1, (arr) + (stack[i - 1].start - 1) * asize, &r)); 754676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 7559566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize)); 7569566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize)); 7579566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 758676f2a66SJacob Faibussowitsch } else { 7599566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize)); 7609566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize)); 7619566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 762676f2a66SJacob Faibussowitsch } 763676f2a66SJacob Faibussowitsch } 764676f2a66SJacob Faibussowitsch /* Update A with merge */ 765676f2a66SJacob Faibussowitsch stack[i - 2].size += stack[i - 1].size; 766676f2a66SJacob Faibussowitsch /* Push C up the stack */ 767676f2a66SJacob Faibussowitsch stack[i - 1].start = stack[i].start; 768676f2a66SJacob Faibussowitsch stack[i - 1].size = stack[i].size; 769676f2a66SJacob Faibussowitsch } else { 770676f2a66SJacob Faibussowitsch /* merge C into B */ 771676f2a66SJacob Faibussowitsch mergeBC: 772676f2a66SJacob Faibussowitsch /* If B[-1] <= C[0] then... you know the drill */ 7734d3610e3SJacob Faibussowitsch if ((*cmp)(arr + (stack[i].start - 1) * asize, arr + (stack[i].start) * asize, ctx) > 0) { 774676f2a66SJacob Faibussowitsch m = stack[i].start; 775676f2a66SJacob Faibussowitsch /* Search B for C[0] insertion */ 7769566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 1].start, stack[i].start - 1, arr + stack[i].start * asize, &l)); 777676f2a66SJacob Faibussowitsch /* Search C for B[-1] insertion */ 7789566063dSJacob Faibussowitsch PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[i].start, stack[i].start + stack[i].size - 1, (arr) + (stack[i].start - 1) * asize, &r)); 779676f2a66SJacob Faibussowitsch if (m - l <= r - m) { 7809566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize)); 7819566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize)); 7829566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 783676f2a66SJacob Faibussowitsch } else { 7849566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize)); 7859566063dSJacob Faibussowitsch PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize)); 7869566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r)); 787676f2a66SJacob Faibussowitsch } 788676f2a66SJacob Faibussowitsch } 789676f2a66SJacob Faibussowitsch /* Update B with merge */ 790676f2a66SJacob Faibussowitsch stack[i - 1].size += stack[i].size; 791676f2a66SJacob Faibussowitsch } 792676f2a66SJacob Faibussowitsch --i; 793676f2a66SJacob Faibussowitsch } else if (stack[i - 1].size <= stack[i].size) { 794676f2a66SJacob Faibussowitsch /* merge C into B */ 795676f2a66SJacob Faibussowitsch goto mergeBC; 796676f2a66SJacob Faibussowitsch } 797676f2a66SJacob Faibussowitsch } 798676f2a66SJacob Faibussowitsch if (itemp == i) break; 799676f2a66SJacob Faibussowitsch } 800676f2a66SJacob Faibussowitsch *stacksize = i; 8013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 802676f2a66SJacob Faibussowitsch } 803676f2a66SJacob Faibussowitsch 804676f2a66SJacob Faibussowitsch /* March sequentially through the array building up a "run" of weakly increasing or strictly decreasing contiguous 805676f2a66SJacob Faibussowitsch elements. Decreasing runs are reversed by swapping. If the run is less than minrun, artificially extend it via either 806676f2a66SJacob Faibussowitsch binary insertion sort or regulat insertion sort */ 807d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortBuildRun_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend) 808d71ae5a4SJacob Faibussowitsch { 809676f2a66SJacob Faibussowitsch const PetscInt re = PetscMin(runstart + minrun, n - 1); 810676f2a66SJacob Faibussowitsch PetscInt ri = runstart; 811676f2a66SJacob Faibussowitsch 812676f2a66SJacob Faibussowitsch PetscFunctionBegin; 8139371c9d4SSatish Balay if (PetscUnlikely(runstart == n - 1)) { 8149371c9d4SSatish Balay *runend = runstart; 8153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8169371c9d4SSatish Balay } 817676f2a66SJacob Faibussowitsch /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */ 8184d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) < 0) { 819676f2a66SJacob Faibussowitsch ++ri; 820676f2a66SJacob Faibussowitsch while (ri < n - 1) { 8214d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) >= 0) break; 822676f2a66SJacob Faibussowitsch ++ri; 823676f2a66SJacob Faibussowitsch } 824676f2a66SJacob Faibussowitsch { 825676f2a66SJacob Faibussowitsch PetscInt lo = runstart, hi = ri; 826ad540459SPierre Jolivet for (; lo < hi; ++lo, --hi) COPYSWAPPY(arr + lo * size, arr + hi * size, tarr, size); 827676f2a66SJacob Faibussowitsch } 828676f2a66SJacob Faibussowitsch } else { 829676f2a66SJacob Faibussowitsch ++ri; 830676f2a66SJacob Faibussowitsch while (ri < n - 1) { 8314d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) < 0) break; 832676f2a66SJacob Faibussowitsch ++ri; 833676f2a66SJacob Faibussowitsch } 834676f2a66SJacob Faibussowitsch } 835676f2a66SJacob Faibussowitsch if (ri < re) { 836676f2a66SJacob Faibussowitsch /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try 837676f2a66SJacob Faibussowitsch binary search */ 838676f2a66SJacob Faibussowitsch if (ri - runstart <= minrun >> 1) { 839676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */ 8403ba16761SJacob Faibussowitsch PetscCall(PetscInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re)); 841676f2a66SJacob Faibussowitsch } else { 8423ba16761SJacob Faibussowitsch PetscCall(PetscBinaryInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re)); 843676f2a66SJacob Faibussowitsch } 844676f2a66SJacob Faibussowitsch *runend = re; 845676f2a66SJacob Faibussowitsch } else *runend = ri; 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 847676f2a66SJacob Faibussowitsch } 848676f2a66SJacob Faibussowitsch 849d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortBuildRunWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend) 850d71ae5a4SJacob Faibussowitsch { 851676f2a66SJacob Faibussowitsch const PetscInt re = PetscMin(runstart + minrun, n - 1); 852676f2a66SJacob Faibussowitsch PetscInt ri = runstart; 853676f2a66SJacob Faibussowitsch 854676f2a66SJacob Faibussowitsch PetscFunctionBegin; 8559371c9d4SSatish Balay if (PetscUnlikely(runstart == n - 1)) { 8569371c9d4SSatish Balay *runend = runstart; 8573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8589371c9d4SSatish Balay } 859676f2a66SJacob Faibussowitsch /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */ 8604d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * asize, arr + (ri * asize), ctx) < 0) { 861676f2a66SJacob Faibussowitsch ++ri; 862676f2a66SJacob Faibussowitsch while (ri < n - 1) { 8634d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * asize, (arr) + ri * asize, ctx) >= 0) break; 864676f2a66SJacob Faibussowitsch ++ri; 865676f2a66SJacob Faibussowitsch } 866676f2a66SJacob Faibussowitsch { 867676f2a66SJacob Faibussowitsch PetscInt lo = runstart, hi = ri; 868ad540459SPierre Jolivet for (; lo < hi; ++lo, --hi) COPYSWAPPY2(arr + lo * asize, arr + hi * asize, asize, barr + lo * bsize, barr + hi * bsize, bsize, atarr); 869676f2a66SJacob Faibussowitsch } 870676f2a66SJacob Faibussowitsch } else { 871676f2a66SJacob Faibussowitsch ++ri; 872676f2a66SJacob Faibussowitsch while (ri < n - 1) { 8734d3610e3SJacob Faibussowitsch if ((*cmp)((arr) + (ri + 1) * asize, (arr) + ri * asize, ctx) < 0) break; 874676f2a66SJacob Faibussowitsch ++ri; 875676f2a66SJacob Faibussowitsch } 876676f2a66SJacob Faibussowitsch } 877676f2a66SJacob Faibussowitsch if (ri < re) { 878676f2a66SJacob Faibussowitsch /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try 879676f2a66SJacob Faibussowitsch binary search */ 880676f2a66SJacob Faibussowitsch if (ri - runstart <= minrun >> 1) { 881676f2a66SJacob Faibussowitsch ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */ 8823ba16761SJacob Faibussowitsch PetscCall(PetscInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re)); 883676f2a66SJacob Faibussowitsch } else { 8843ba16761SJacob Faibussowitsch PetscCall(PetscBinaryInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re)); 885676f2a66SJacob Faibussowitsch } 886676f2a66SJacob Faibussowitsch *runend = re; 887676f2a66SJacob Faibussowitsch } else *runend = ri; 8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 889676f2a66SJacob Faibussowitsch } 890676f2a66SJacob Faibussowitsch 8914d3610e3SJacob Faibussowitsch /*@C 8921d27aa22SBarry Smith PetscTimSort - Sorts an array in place in increasing order using Tim Peters <https://bugs.python.org/file4451/timsort.txt> adaptive sorting algorithm. 893676f2a66SJacob Faibussowitsch 894cc4c1da9SBarry Smith Not Collective, No Fortran Support 895676f2a66SJacob Faibussowitsch 896676f2a66SJacob Faibussowitsch Input Parameters: 897676f2a66SJacob Faibussowitsch + n - number of values 898676f2a66SJacob Faibussowitsch . arr - array to be sorted 899676f2a66SJacob Faibussowitsch . size - size in bytes of the datatype held in arr 9004d3610e3SJacob Faibussowitsch . cmp - function pointer to comparison function 9014d3610e3SJacob Faibussowitsch - ctx - optional context to be passed to comparison function, NULL if not needed 902676f2a66SJacob Faibussowitsch 9032fe279fdSBarry Smith Output Parameter: 904676f2a66SJacob Faibussowitsch . arr - sorted array 905676f2a66SJacob Faibussowitsch 9061d27aa22SBarry Smith Level: developer 9071d27aa22SBarry Smith 9088da36234SJacob Faibussowitsch Notes: 9098da36234SJacob Faibussowitsch Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous 9108da36234SJacob Faibussowitsch sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims to 9118da36234SJacob Faibussowitsch select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced arrays. To 9128da36234SJacob Faibussowitsch do so it repeatedly triggers attempts throughout to merge adjacent runs. 9138da36234SJacob Faibussowitsch 9148da36234SJacob Faibussowitsch Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search 9158da36234SJacob Faibussowitsch the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in 9168da36234SJacob Faibussowitsch bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these 9178da36234SJacob Faibussowitsch searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they 9188da36234SJacob Faibussowitsch likely all contain similar data. 9198da36234SJacob Faibussowitsch 920aec76313SJacob Faibussowitsch Example Usage: 921811af0c4SBarry Smith The comparison function must follow the `qsort()` comparison function paradigm, returning the sign of the difference 922676f2a66SJacob Faibussowitsch between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user 923676f2a66SJacob Faibussowitsch may also 924676f2a66SJacob Faibussowitsch change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if 925676f2a66SJacob Faibussowitsch the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing 9268da36234SJacob Faibussowitsch order 927676f2a66SJacob Faibussowitsch 928676f2a66SJacob Faibussowitsch .vb 9291d27aa22SBarry Smith int increasing_comparison_function(const void *left, const void *right, void *ctx) { 930676f2a66SJacob Faibussowitsch my_type l = *(my_type *) left, r = *(my_type *) right; 931676f2a66SJacob Faibussowitsch return (l < r) ? -1 : (l > r); 932676f2a66SJacob Faibussowitsch } 933676f2a66SJacob Faibussowitsch .ve 934aec76313SJacob Faibussowitsch 9354d3610e3SJacob Faibussowitsch Note the context is unused here but you may use it to pass and subsequently access whatever information required 9364d3610e3SJacob Faibussowitsch inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function. 937676f2a66SJacob Faibussowitsch Then pass the function 938676f2a66SJacob Faibussowitsch .vb 9391d27aa22SBarry Smith PetscTimSort(n, arr, sizeof(arr[0]), increasing_comparison_function, ctx) 940676f2a66SJacob Faibussowitsch .ve 941676f2a66SJacob Faibussowitsch 9424d3610e3SJacob Faibussowitsch Fortran Notes: 943*feaf08eaSBarry Smith To use this from Fortran you must write a comparison subroutine with 4 arguments which accepts `left`, `right`, `ctx` and 944*feaf08eaSBarry Smith returns `result`. For example 9454d3610e3SJacob Faibussowitsch .vb 9464d3610e3SJacob Faibussowitsch subroutine CompareIntegers(left,right,ctx,result) 9474d3610e3SJacob Faibussowitsch implicit none 9484d3610e3SJacob Faibussowitsch 9494d3610e3SJacob Faibussowitsch PetscInt,intent(in) :: left, right 9504d3610e3SJacob Faibussowitsch type(UserCtx) :: ctx 9514d3610e3SJacob Faibussowitsch integer,intent(out) :: result 9524d3610e3SJacob Faibussowitsch 9534d3610e3SJacob Faibussowitsch if (left < right) then 9544d3610e3SJacob Faibussowitsch result = -1 9554d3610e3SJacob Faibussowitsch else if (left == right) then 9564d3610e3SJacob Faibussowitsch result = 0 9574d3610e3SJacob Faibussowitsch else 9584d3610e3SJacob Faibussowitsch result = 1 9594d3610e3SJacob Faibussowitsch end if 9604d3610e3SJacob Faibussowitsch return 9614d3610e3SJacob Faibussowitsch end subroutine CompareIntegers 9624d3610e3SJacob Faibussowitsch .ve 9634d3610e3SJacob Faibussowitsch 964db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscIntSortSemiOrdered()`, `PetscRealSortSemiOrdered()`, `PetscMPIIntSortSemiOrdered()` 9654d3610e3SJacob Faibussowitsch @*/ 966d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *, void *), void *ctx) 967d71ae5a4SJacob Faibussowitsch { 968676f2a66SJacob Faibussowitsch PetscInt stacksize = 0, minrun, runstart = 0, runend = 0; 969676f2a66SJacob Faibussowitsch PetscTimSortStack runstack[128]; 970676f2a66SJacob Faibussowitsch PetscTimSortBuffer buff; 971676f2a66SJacob Faibussowitsch /* stacksize = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements. 972676f2a66SJacob Faibussowitsch It is so unlikely that this limit is reached that this is __never__ checked for */ 973676f2a66SJacob Faibussowitsch 974676f2a66SJacob Faibussowitsch PetscFunctionBegin; 975676f2a66SJacob Faibussowitsch /* Compute minrun. Minrun should be (32, 65) such that N/minrun 976676f2a66SJacob Faibussowitsch is a power of 2 or one plus a power of 2 */ 977676f2a66SJacob Faibussowitsch { 978676f2a66SJacob Faibussowitsch PetscInt t = n, r = 0; 979676f2a66SJacob Faibussowitsch /* r becomes 1 if the least significant bits contain at least one off bit */ 980676f2a66SJacob Faibussowitsch while (t >= 64) { 981676f2a66SJacob Faibussowitsch r |= t & 1; 982676f2a66SJacob Faibussowitsch t >>= 1; 983676f2a66SJacob Faibussowitsch } 984676f2a66SJacob Faibussowitsch minrun = t + r; 985676f2a66SJacob Faibussowitsch } 9864d3610e3SJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) { 9879566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun)); 9884d3610e3SJacob Faibussowitsch if (n < 64) { 9899566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "n %" PetscInt_FMT " < 64, consider using PetscSortInt() instead\n", n)); 99008401ef6SPierre Jolivet } else PetscCheck(!(minrun < 32) && !(minrun > 65), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Calculated minrun %" PetscInt_FMT " not in range (32,65)", minrun); 9914d3610e3SJacob Faibussowitsch } 9929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((size_t)minrun * size, &buff.ptr)); 993676f2a66SJacob Faibussowitsch buff.size = (size_t)minrun * size; 994676f2a66SJacob Faibussowitsch buff.maxsize = (size_t)n * size; 995676f2a66SJacob Faibussowitsch MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL; 996676f2a66SJacob Faibussowitsch while (runstart < n) { 997676f2a66SJacob Faibussowitsch /* Check if additional entries are at least partially ordered and build natural run */ 9989566063dSJacob Faibussowitsch PetscCall(PetscTimSortBuildRun_Private((char *)arr, buff.ptr, size, cmp, ctx, n, minrun, runstart, &runend)); 999676f2a66SJacob Faibussowitsch runstack[stacksize].start = runstart; 1000676f2a66SJacob Faibussowitsch runstack[stacksize].size = runend - runstart + 1; 10019566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, &stacksize)); 1002676f2a66SJacob Faibussowitsch ++stacksize; 1003676f2a66SJacob Faibussowitsch runstart = runend + 1; 1004676f2a66SJacob Faibussowitsch } 1005676f2a66SJacob Faibussowitsch /* Have been inside while, so discard last stacksize++ */ 1006676f2a66SJacob Faibussowitsch --stacksize; 10079566063dSJacob Faibussowitsch PetscCall(PetscTimSortForceCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, stacksize)); 10089566063dSJacob Faibussowitsch PetscCall(PetscFree(buff.ptr)); 1009676f2a66SJacob Faibussowitsch MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL; 10103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1011676f2a66SJacob Faibussowitsch } 1012676f2a66SJacob Faibussowitsch 10134d3610e3SJacob Faibussowitsch /*@C 10141d27aa22SBarry Smith PetscTimSortWithArray - Sorts an array in place in increasing order using Tim Peters <https://bugs.python.org/file4451/timsort.txt> adaptive sorting algorithm and 1015676f2a66SJacob Faibussowitsch reorders a second array to match the first. The arrays need not be the same type. 1016676f2a66SJacob Faibussowitsch 1017cc4c1da9SBarry Smith Not Collective, No Fortran Support 1018676f2a66SJacob Faibussowitsch 1019676f2a66SJacob Faibussowitsch Input Parameters: 1020676f2a66SJacob Faibussowitsch + n - number of values 1021676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in arr 10226b867d5aSJose E. Roman . bsize - size in bytes of the datatype held in barr 10234d3610e3SJacob Faibussowitsch . cmp - function pointer to comparison function 10244d3610e3SJacob Faibussowitsch - ctx - optional context to be passed to comparison function, NULL if not needed 1025676f2a66SJacob Faibussowitsch 10266b867d5aSJose E. Roman Input/Output Parameters: 10276b867d5aSJose E. Roman + arr - array to be sorted, on output it is sorted 10286b867d5aSJose E. Roman - barr - array to be reordered, on output it is reordered 10298da36234SJacob Faibussowitsch 10301d27aa22SBarry Smith Level: developer 10311d27aa22SBarry Smith 10328da36234SJacob Faibussowitsch Notes: 10338da36234SJacob Faibussowitsch The arrays need not be of the same type, however barr MUST contain at least as many elements as arr and the two CANNOT 10348da36234SJacob Faibussowitsch overlap. 10358da36234SJacob Faibussowitsch 10368da36234SJacob Faibussowitsch Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous 10378da36234SJacob Faibussowitsch sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims 10388da36234SJacob Faibussowitsch to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced 10398da36234SJacob Faibussowitsch arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs. 10408da36234SJacob Faibussowitsch 10418da36234SJacob Faibussowitsch Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search 10428da36234SJacob Faibussowitsch the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in 10438da36234SJacob Faibussowitsch bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these 10448da36234SJacob Faibussowitsch searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they 10458da36234SJacob Faibussowitsch likely all contain similar data. 1046676f2a66SJacob Faibussowitsch 1047aec76313SJacob Faibussowitsch Example Usage: 1048811af0c4SBarry Smith The comparison function must follow the `qsort()` comparison function paradigm, returning the sign of the difference 10491d27aa22SBarry Smith between its arguments. If left < right \: return -1, if left == right \: return 0, if left > right \: return 1. The user 1050676f2a66SJacob Faibussowitsch may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only 1051676f2a66SJacob Faibussowitsch guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in 10528da36234SJacob Faibussowitsch increasing order 1053676f2a66SJacob Faibussowitsch 1054676f2a66SJacob Faibussowitsch .vb 10551d27aa22SBarry Smith int increasing_comparison_function(const void *left, const void *right, void *ctx) { 1056676f2a66SJacob Faibussowitsch my_type l = *(my_type *) left, r = *(my_type *) right; 1057676f2a66SJacob Faibussowitsch return (l < r) ? -1 : (l > r); 1058676f2a66SJacob Faibussowitsch } 1059676f2a66SJacob Faibussowitsch .ve 1060aec76313SJacob Faibussowitsch 10614d3610e3SJacob Faibussowitsch Note the context is unused here but you may use it to pass and subsequently access whatever information required 10624d3610e3SJacob Faibussowitsch inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function. 1063676f2a66SJacob Faibussowitsch Then pass the function 1064676f2a66SJacob Faibussowitsch .vb 10651d27aa22SBarry Smith PetscTimSortWithArray(n, arr, sizeof(arr[0]), barr, sizeof(barr[0]), increasing_comparison_function, ctx) 1066676f2a66SJacob Faibussowitsch .ve 1067676f2a66SJacob Faibussowitsch 10684d3610e3SJacob Faibussowitsch Fortran Notes: 1069*feaf08eaSBarry Smith To use this from Fortran you must write a comparison subroutine with 4 arguments which accepts `left`, `right`, `ctx` and 1070*feaf08eaSBarry Smith returns `result`. For example 10714d3610e3SJacob Faibussowitsch .vb 10724d3610e3SJacob Faibussowitsch subroutine CompareIntegers(left,right,ctx,result) 10734d3610e3SJacob Faibussowitsch implicit none 10744d3610e3SJacob Faibussowitsch 10754d3610e3SJacob Faibussowitsch PetscInt,intent(in) :: left, right 10764d3610e3SJacob Faibussowitsch type(UserCtx) :: ctx 10774d3610e3SJacob Faibussowitsch integer,intent(out) :: result 10784d3610e3SJacob Faibussowitsch 10794d3610e3SJacob Faibussowitsch if (left < right) then 10804d3610e3SJacob Faibussowitsch result = -1 10814d3610e3SJacob Faibussowitsch else if (left == right) then 10824d3610e3SJacob Faibussowitsch result = 0 10834d3610e3SJacob Faibussowitsch else 10844d3610e3SJacob Faibussowitsch result = 1 10854d3610e3SJacob Faibussowitsch end if 10864d3610e3SJacob Faibussowitsch return 10874d3610e3SJacob Faibussowitsch end subroutine CompareIntegers 10884d3610e3SJacob Faibussowitsch .ve 10894d3610e3SJacob Faibussowitsch 1090db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscIntSortSemiOrderedWithArray()`, `PetscRealSortSemiOrderedWithArrayInt()`, `PetscMPIIntSortSemiOrderedWithArray()` 10914d3610e3SJacob Faibussowitsch @*/ 1092d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscTimSortWithArray(PetscInt n, void *arr, size_t asize, void *barr, size_t bsize, int (*cmp)(const void *, const void *, void *), void *ctx) 1093d71ae5a4SJacob Faibussowitsch { 1094676f2a66SJacob Faibussowitsch PetscInt stacksize = 0, minrun, runstart = 0, runend = 0; 1095676f2a66SJacob Faibussowitsch PetscTimSortStack runstack[128]; 1096676f2a66SJacob Faibussowitsch PetscTimSortBuffer abuff, bbuff; 1097676f2a66SJacob Faibussowitsch /* stacksize = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements. 1098676f2a66SJacob Faibussowitsch It is so unlikely that this limit is reached that this is __never__ checked for */ 1099676f2a66SJacob Faibussowitsch 1100676f2a66SJacob Faibussowitsch PetscFunctionBegin; 1101676f2a66SJacob Faibussowitsch /* Compute minrun. Minrun should be (32, 65) such that N/minrun 1102676f2a66SJacob Faibussowitsch is a power of 2 or one plus a power of 2 */ 1103676f2a66SJacob Faibussowitsch { 1104676f2a66SJacob Faibussowitsch PetscInt t = n, r = 0; 1105676f2a66SJacob Faibussowitsch /* r becomes 1 if the least significant bits contain at least one off bit */ 1106676f2a66SJacob Faibussowitsch while (t >= 64) { 1107676f2a66SJacob Faibussowitsch r |= t & 1; 1108676f2a66SJacob Faibussowitsch t >>= 1; 1109676f2a66SJacob Faibussowitsch } 1110676f2a66SJacob Faibussowitsch minrun = t + r; 1111676f2a66SJacob Faibussowitsch } 11124d3610e3SJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) { 11139566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun)); 1114cc73adaaSBarry Smith PetscCheck(n < 64 || (minrun >= 32 && minrun <= 65), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Calculated minrun %" PetscInt_FMT " not in range (32,65)", minrun); 11154d3610e3SJacob Faibussowitsch } 11169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((size_t)minrun * asize, &abuff.ptr)); 1117676f2a66SJacob Faibussowitsch abuff.size = (size_t)minrun * asize; 1118676f2a66SJacob Faibussowitsch abuff.maxsize = (size_t)n * asize; 11199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((size_t)minrun * bsize, &bbuff.ptr)); 1120676f2a66SJacob Faibussowitsch bbuff.size = (size_t)minrun * bsize; 1121676f2a66SJacob Faibussowitsch bbuff.maxsize = (size_t)n * bsize; 1122676f2a66SJacob Faibussowitsch MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL; 1123676f2a66SJacob Faibussowitsch while (runstart < n) { 1124676f2a66SJacob Faibussowitsch /* Check if additional entries are at least partially ordered and build natural run */ 11259566063dSJacob Faibussowitsch PetscCall(PetscTimSortBuildRunWithArray_Private((char *)arr, abuff.ptr, asize, (char *)barr, bbuff.ptr, bsize, cmp, ctx, n, minrun, runstart, &runend)); 1126676f2a66SJacob Faibussowitsch runstack[stacksize].start = runstart; 1127676f2a66SJacob Faibussowitsch runstack[stacksize].size = runend - runstart + 1; 11289566063dSJacob Faibussowitsch PetscCall(PetscTimSortMergeCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, &stacksize)); 1129676f2a66SJacob Faibussowitsch ++stacksize; 1130676f2a66SJacob Faibussowitsch runstart = runend + 1; 1131676f2a66SJacob Faibussowitsch } 1132676f2a66SJacob Faibussowitsch /* Have been inside while, so discard last stacksize++ */ 1133676f2a66SJacob Faibussowitsch --stacksize; 11349566063dSJacob Faibussowitsch PetscCall(PetscTimSortForceCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, stacksize)); 11359566063dSJacob Faibussowitsch PetscCall(PetscFree(abuff.ptr)); 11369566063dSJacob Faibussowitsch PetscCall(PetscFree(bbuff.ptr)); 1137676f2a66SJacob Faibussowitsch MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL; 11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1139676f2a66SJacob Faibussowitsch } 1140676f2a66SJacob Faibussowitsch 1141676f2a66SJacob Faibussowitsch /*@ 1142811af0c4SBarry Smith PetscIntSortSemiOrdered - Sorts an array of `PetscInt` in place in increasing order. 1143676f2a66SJacob Faibussowitsch 1144676f2a66SJacob Faibussowitsch Not Collective 1145676f2a66SJacob Faibussowitsch 1146676f2a66SJacob Faibussowitsch Input Parameters: 1147676f2a66SJacob Faibussowitsch + n - number of values 1148676f2a66SJacob Faibussowitsch - arr - array of integers 1149676f2a66SJacob Faibussowitsch 11502fe279fdSBarry Smith Output Parameter: 1151676f2a66SJacob Faibussowitsch . arr - sorted array of integers 1152676f2a66SJacob Faibussowitsch 11531d27aa22SBarry Smith Level: intermediate 11541d27aa22SBarry Smith 1155676f2a66SJacob Faibussowitsch Notes: 1156811af0c4SBarry Smith If the array is less than 64 entries long `PetscSortInt()` is automatically used. 1157676f2a66SJacob Faibussowitsch 1158811af0c4SBarry Smith This function serves as an alternative to `PetscSortInt()`. While this function works for any array of integers it is 1159676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1160a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1161676f2a66SJacob Faibussowitsch 1162db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortInt()`, `PetscSortIntWithPermutation()` 1163676f2a66SJacob Faibussowitsch @*/ 1164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrdered(PetscInt n, PetscInt arr[]) 1165d71ae5a4SJacob Faibussowitsch { 1166676f2a66SJacob Faibussowitsch PetscFunctionBegin; 11673ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 11684f572ea9SToby Isaac PetscAssertPointer(arr, 2); 1169676f2a66SJacob Faibussowitsch if (n < 64) { 11709566063dSJacob Faibussowitsch PetscCall(PetscSortInt(n, arr)); 1171676f2a66SJacob Faibussowitsch } else { 11729566063dSJacob Faibussowitsch PetscCall(PetscTimSort(n, arr, sizeof(PetscInt), Compare_PetscInt_Private, NULL)); 1173676f2a66SJacob Faibussowitsch } 11743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1175676f2a66SJacob Faibussowitsch } 1176676f2a66SJacob Faibussowitsch 1177676f2a66SJacob Faibussowitsch /*@ 1178811af0c4SBarry Smith PetscIntSortSemiOrderedWithArray - Sorts an array of `PetscInt` in place in increasing order and reorders a second 1179811af0c4SBarry Smith `PetscInt` array to match the first. 1180676f2a66SJacob Faibussowitsch 1181676f2a66SJacob Faibussowitsch Not Collective 1182676f2a66SJacob Faibussowitsch 11836b867d5aSJose E. Roman Input Parameter: 11846b867d5aSJose E. Roman . n - number of values 1185676f2a66SJacob Faibussowitsch 11866b867d5aSJose E. Roman Input/Output Parameters: 11876b867d5aSJose E. Roman + arr1 - array of integers to be sorted, modified on output 11886b867d5aSJose E. Roman - arr2 - array of integers to be reordered, modified on output 1189676f2a66SJacob Faibussowitsch 11901d27aa22SBarry Smith Level: intermediate 11911d27aa22SBarry Smith 1192676f2a66SJacob Faibussowitsch Notes: 1193676f2a66SJacob Faibussowitsch The arrays CANNOT overlap. 1194676f2a66SJacob Faibussowitsch 11951d27aa22SBarry Smith This function serves as an alternative to `PetscSortIntWithArray()`. While this function works for any array of integers it is 1196676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1197a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1198676f2a66SJacob Faibussowitsch 1199db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortIntWithArray()`, `PetscSortIntWithPermutation()` 1200676f2a66SJacob Faibussowitsch @*/ 1201d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt n, PetscInt arr1[], PetscInt arr2[]) 1202d71ae5a4SJacob Faibussowitsch { 1203676f2a66SJacob Faibussowitsch PetscFunctionBegin; 12043ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 12054f572ea9SToby Isaac PetscAssertPointer(arr1, 2); 12064f572ea9SToby Isaac PetscAssertPointer(arr2, 3); 12076e260675SJacob Faibussowitsch /* cannot export out to PetscIntSortWithArray here since it isn't stable */ 12089566063dSJacob Faibussowitsch PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscInt), arr2, sizeof(PetscInt), Compare_PetscInt_Private, NULL)); 12093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1210676f2a66SJacob Faibussowitsch } 1211676f2a66SJacob Faibussowitsch 1212676f2a66SJacob Faibussowitsch /*@ 1213811af0c4SBarry Smith PetscMPIIntSortSemiOrdered - Sorts an array of `PetscMPIInt` in place in increasing order. 1214676f2a66SJacob Faibussowitsch 1215676f2a66SJacob Faibussowitsch Not Collective 1216676f2a66SJacob Faibussowitsch 1217676f2a66SJacob Faibussowitsch Input Parameters: 1218676f2a66SJacob Faibussowitsch + n - number of values 1219811af0c4SBarry Smith - arr - array of `PetscMPIInt` 1220676f2a66SJacob Faibussowitsch 12212fe279fdSBarry Smith Output Parameter: 1222676f2a66SJacob Faibussowitsch . arr - sorted array of integers 1223676f2a66SJacob Faibussowitsch 12241d27aa22SBarry Smith Level: intermediate 12251d27aa22SBarry Smith 1226676f2a66SJacob Faibussowitsch Notes: 1227811af0c4SBarry Smith If the array is less than 64 entries long `PetscSortMPIInt()` is automatically used. 1228676f2a66SJacob Faibussowitsch 1229811af0c4SBarry Smith This function serves as an alternative to `PetscSortMPIInt()`. While this function works for any array of `PetscMPIInt` it is 1230676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1231a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1232676f2a66SJacob Faibussowitsch 1233db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortMPIInt()` 1234676f2a66SJacob Faibussowitsch @*/ 1235d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt n, PetscMPIInt arr[]) 1236d71ae5a4SJacob Faibussowitsch { 1237676f2a66SJacob Faibussowitsch PetscFunctionBegin; 12383ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 12394f572ea9SToby Isaac PetscAssertPointer(arr, 2); 1240676f2a66SJacob Faibussowitsch if (n < 64) { 12419566063dSJacob Faibussowitsch PetscCall(PetscSortMPIInt(n, arr)); 1242676f2a66SJacob Faibussowitsch } else { 12439566063dSJacob Faibussowitsch PetscCall(PetscTimSort(n, arr, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL)); 1244676f2a66SJacob Faibussowitsch } 12453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1246676f2a66SJacob Faibussowitsch } 1247676f2a66SJacob Faibussowitsch 1248676f2a66SJacob Faibussowitsch /*@ 1249811af0c4SBarry Smith PetscMPIIntSortSemiOrderedWithArray - Sorts an array of `PetscMPIInt` in place in increasing order and reorders a second `PetscMPIInt` 1250676f2a66SJacob Faibussowitsch array to match the first. 1251676f2a66SJacob Faibussowitsch 1252676f2a66SJacob Faibussowitsch Not Collective 1253676f2a66SJacob Faibussowitsch 12546b867d5aSJose E. Roman Input Parameter: 12556b867d5aSJose E. Roman . n - number of values 1256676f2a66SJacob Faibussowitsch 12576b867d5aSJose E. Roman Input/Output Parameters: 1258aec76313SJacob Faibussowitsch + arr1 - array of integers to be sorted, modified on output 12596b867d5aSJose E. Roman - arr2 - array of integers to be reordered, modified on output 1260676f2a66SJacob Faibussowitsch 12611d27aa22SBarry Smith Level: intermediate 12621d27aa22SBarry Smith 1263676f2a66SJacob Faibussowitsch Notes: 1264676f2a66SJacob Faibussowitsch The arrays CANNOT overlap. 1265676f2a66SJacob Faibussowitsch 1266811af0c4SBarry Smith This function serves as an alternative to `PetscSortMPIIntWithArray()`. While this function works for any array of integers it is 1267676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1268a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1269676f2a66SJacob Faibussowitsch 1270db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortMPIIntWithArray()`, `PetscSortMPIIntWithPermutation()` 1271676f2a66SJacob Faibussowitsch @*/ 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt n, PetscMPIInt arr1[], PetscMPIInt arr2[]) 1273d71ae5a4SJacob Faibussowitsch { 1274676f2a66SJacob Faibussowitsch PetscFunctionBegin; 12753ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 12764f572ea9SToby Isaac PetscAssertPointer(arr1, 2); 12774f572ea9SToby Isaac PetscAssertPointer(arr2, 3); 12786e260675SJacob Faibussowitsch /* cannot export out to PetscMPIIntSortWithArray here since it isn't stable */ 12799566063dSJacob Faibussowitsch PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscMPIInt), arr2, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL)); 12803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1281676f2a66SJacob Faibussowitsch } 1282676f2a66SJacob Faibussowitsch 1283676f2a66SJacob Faibussowitsch /*@ 1284811af0c4SBarry Smith PetscRealSortSemiOrdered - Sorts an array of `PetscReal` in place in increasing order. 1285676f2a66SJacob Faibussowitsch 1286676f2a66SJacob Faibussowitsch Not Collective 1287676f2a66SJacob Faibussowitsch 1288676f2a66SJacob Faibussowitsch Input Parameters: 1289676f2a66SJacob Faibussowitsch + n - number of values 1290811af0c4SBarry Smith - arr - array of `PetscReal` 1291676f2a66SJacob Faibussowitsch 12922fe279fdSBarry Smith Output Parameter: 1293676f2a66SJacob Faibussowitsch . arr - sorted array of integers 1294676f2a66SJacob Faibussowitsch 12951d27aa22SBarry Smith Level: intermediate 12961d27aa22SBarry Smith 1297676f2a66SJacob Faibussowitsch Notes: 1298811af0c4SBarry Smith If the array is less than 64 entries long `PetscSortReal()` is automatically used. 1299676f2a66SJacob Faibussowitsch 1300811af0c4SBarry Smith This function serves as an alternative to `PetscSortReal()`. While this function works for any array of `PetscReal` it is 1301676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1302a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1303676f2a66SJacob Faibussowitsch 1304db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortReal()`, `PetscSortRealWithPermutation()` 1305676f2a66SJacob Faibussowitsch @*/ 1306d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrdered(PetscInt n, PetscReal arr[]) 1307d71ae5a4SJacob Faibussowitsch { 1308676f2a66SJacob Faibussowitsch PetscFunctionBegin; 13093ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 13104f572ea9SToby Isaac PetscAssertPointer(arr, 2); 1311676f2a66SJacob Faibussowitsch if (n < 64) { 13129566063dSJacob Faibussowitsch PetscCall(PetscSortReal(n, arr)); 1313676f2a66SJacob Faibussowitsch } else { 13149566063dSJacob Faibussowitsch PetscCall(PetscTimSort(n, arr, sizeof(PetscReal), Compare_PetscReal_Private, NULL)); 1315676f2a66SJacob Faibussowitsch } 13163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1317676f2a66SJacob Faibussowitsch } 1318676f2a66SJacob Faibussowitsch 1319676f2a66SJacob Faibussowitsch /*@ 1320811af0c4SBarry Smith PetscRealSortSemiOrderedWithArrayInt - Sorts an array of `PetscReal` in place in increasing order and reorders a second 1321811af0c4SBarry Smith array of `PetscInt` to match the first. 1322676f2a66SJacob Faibussowitsch 1323676f2a66SJacob Faibussowitsch Not Collective 1324676f2a66SJacob Faibussowitsch 13256b867d5aSJose E. Roman Input Parameter: 13266b867d5aSJose E. Roman . n - number of values 1327676f2a66SJacob Faibussowitsch 13286b867d5aSJose E. Roman Input/Output Parameters: 1329aec76313SJacob Faibussowitsch + arr1 - array of `PetscReal` to be sorted, modified on output 1330811af0c4SBarry Smith - arr2 - array of `PetscInt` to be reordered, modified on output 1331676f2a66SJacob Faibussowitsch 13321d27aa22SBarry Smith Level: intermediate 13331d27aa22SBarry Smith 1334676f2a66SJacob Faibussowitsch Notes: 1335811af0c4SBarry Smith This function serves as an alternative to `PetscSortRealWithArray()`. While this function works for any array of `PetscReal` it is 1336676f2a66SJacob Faibussowitsch significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__ 1337a5b23f4aSJose E. Roman recommended that the user benchmark their code to see which routine is fastest. 1338676f2a66SJacob Faibussowitsch 1339db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortRealWithArrayInt()`, `PetscSortRealWithPermutation()` 1340676f2a66SJacob Faibussowitsch @*/ 1341d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt n, PetscReal arr1[], PetscInt arr2[]) 1342d71ae5a4SJacob Faibussowitsch { 1343676f2a66SJacob Faibussowitsch PetscFunctionBegin; 13443ba16761SJacob Faibussowitsch if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS); 13454f572ea9SToby Isaac PetscAssertPointer(arr1, 2); 13464f572ea9SToby Isaac PetscAssertPointer(arr2, 3); 13476e260675SJacob Faibussowitsch /* cannot export out to PetscRealSortWithArrayInt here since it isn't stable */ 13489566063dSJacob Faibussowitsch PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscReal), arr2, sizeof(PetscInt), Compare_PetscReal_Private, NULL)); 13493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1350676f2a66SJacob Faibussowitsch } 1351