xref: /petsc/src/sys/utils/sortso.c (revision 811af0c4b09a35de4306c442f88bd09fdc09897d)
1676f2a66SJacob Faibussowitsch #include <petscsys.h> /*I  "petscsys.h"  I*/
2676f2a66SJacob Faibussowitsch #include <petsc/private/petscimpl.h>
3676f2a66SJacob Faibussowitsch 
49371c9d4SSatish Balay static inline int Compare_PetscMPIInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) {
5676f2a66SJacob Faibussowitsch   PetscMPIInt l = *(PetscMPIInt *)left, r = *(PetscMPIInt *)right;
6676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
7676f2a66SJacob Faibussowitsch }
8676f2a66SJacob Faibussowitsch 
99371c9d4SSatish Balay static inline int Compare_PetscInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) {
10676f2a66SJacob Faibussowitsch   PetscInt l = *(PetscInt *)left, r = *(PetscInt *)right;
11676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
12676f2a66SJacob Faibussowitsch }
13676f2a66SJacob Faibussowitsch 
149371c9d4SSatish Balay static inline int Compare_PetscReal_Private(const void *left, const void *right, PETSC_UNUSED void *ctx) {
15676f2a66SJacob Faibussowitsch   PetscReal l = *(PetscReal *)left, r = *(PetscReal *)right;
16676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
17676f2a66SJacob Faibussowitsch }
18676f2a66SJacob Faibussowitsch 
19676f2a66SJacob Faibussowitsch #define MIN_GALLOP_CONST_GLOBAL 8
20676f2a66SJacob Faibussowitsch static PetscInt MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
214d3610e3SJacob Faibussowitsch typedef int (*CompFunc)(const void *, const void *, void *);
22676f2a66SJacob Faibussowitsch 
23676f2a66SJacob Faibussowitsch /* Mostly to force clang uses the builtins, but can't hurt to have gcc compilers also do the same */
24676f2a66SJacob Faibussowitsch #if !defined(PETSC_USE_DEBUG) && defined(__GNUC__)
259371c9d4SSatish Balay static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size) {
26676f2a66SJacob Faibussowitsch   __builtin_memcpy(t, b, size);
27676f2a66SJacob Faibussowitsch   __builtin_memmove(b, a, size);
28676f2a66SJacob Faibussowitsch   __builtin_memcpy(a, t, size);
29676f2a66SJacob Faibussowitsch   return;
30676f2a66SJacob Faibussowitsch }
31676f2a66SJacob Faibussowitsch 
329371c9d4SSatish Balay static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t) {
336e260675SJacob Faibussowitsch   __builtin_memcpy(t, ar, asize);
346e260675SJacob Faibussowitsch   __builtin_memmove(ar, al, asize);
356e260675SJacob Faibussowitsch   __builtin_memcpy(al, t, asize);
366e260675SJacob Faibussowitsch   __builtin_memcpy(t, br, bsize);
376e260675SJacob Faibussowitsch   __builtin_memmove(br, bl, bsize);
386e260675SJacob Faibussowitsch   __builtin_memcpy(bl, t, bsize);
39676f2a66SJacob Faibussowitsch   return;
40676f2a66SJacob Faibussowitsch }
41676f2a66SJacob Faibussowitsch 
429371c9d4SSatish Balay static inline void Petsc_memcpy(char *dest, const char *src, size_t size) {
43676f2a66SJacob Faibussowitsch   __builtin_memcpy(dest, src, size);
44676f2a66SJacob Faibussowitsch   return;
45676f2a66SJacob Faibussowitsch }
46676f2a66SJacob Faibussowitsch 
479371c9d4SSatish Balay static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) {
48676f2a66SJacob Faibussowitsch   __builtin_memcpy(adest, asrc, asize);
496e260675SJacob Faibussowitsch   __builtin_memcpy(bdest, bsrc, bsize);
50676f2a66SJacob Faibussowitsch   return;
51676f2a66SJacob Faibussowitsch }
52676f2a66SJacob Faibussowitsch 
539371c9d4SSatish Balay static inline void Petsc_memmove(char *dest, const char *src, size_t size) {
54676f2a66SJacob Faibussowitsch   __builtin_memmove(dest, src, size);
55676f2a66SJacob Faibussowitsch   return;
56676f2a66SJacob Faibussowitsch }
57676f2a66SJacob Faibussowitsch 
589371c9d4SSatish Balay static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) {
59676f2a66SJacob Faibussowitsch   __builtin_memmove(adest, asrc, asize);
60676f2a66SJacob Faibussowitsch   __builtin_memmove(bdest, bsrc, bsize);
61676f2a66SJacob Faibussowitsch   return;
62676f2a66SJacob Faibussowitsch }
63676f2a66SJacob Faibussowitsch #else
649371c9d4SSatish Balay static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size) {
65676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, b, size));
679566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(b, a, size));
689566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(a, t, size));
69676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
70676f2a66SJacob Faibussowitsch }
71676f2a66SJacob Faibussowitsch 
729371c9d4SSatish Balay static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t) {
73676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
749566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, ar, asize));
759566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(ar, al, asize));
769566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(al, t, asize));
779566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(t, br, bsize));
789566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(br, bl, bsize));
799566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(bl, t, bsize));
80676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
81676f2a66SJacob Faibussowitsch }
82676f2a66SJacob Faibussowitsch 
839371c9d4SSatish Balay static inline void Petsc_memcpy(char *dest, const char *src, size_t size) {
84676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(dest, src, size));
86676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
87676f2a66SJacob Faibussowitsch }
88676f2a66SJacob Faibussowitsch 
899371c9d4SSatish Balay static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) {
90676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(adest, asrc, asize));
929566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemcpy(bdest, bsrc, bsize));
93676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
94676f2a66SJacob Faibussowitsch }
95676f2a66SJacob Faibussowitsch 
969371c9d4SSatish Balay static inline void Petsc_memmove(char *dest, const char *src, size_t size) {
97676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(dest, src, size));
99676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
100676f2a66SJacob Faibussowitsch }
101676f2a66SJacob Faibussowitsch 
1029371c9d4SSatish Balay static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize) {
103676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(adest, asrc, asize));
1059566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscMemmove(bdest, bsrc, bsize));
106676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
107676f2a66SJacob Faibussowitsch }
108676f2a66SJacob Faibussowitsch #endif
109676f2a66SJacob Faibussowitsch 
110676f2a66SJacob 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] >
111676f2a66SJacob Faibussowitsch  x. Output also inclusive.
112676f2a66SJacob Faibussowitsch 
113676f2a66SJacob Faibussowitsch  NOTE: Both gallopsearch functions CANNOT distinguish between inserting AFTER the entire array vs inserting at entry n!! For example for an array:
114676f2a66SJacob Faibussowitsch 
115676f2a66SJacob Faibussowitsch    {0,1,2,3,4,5}
116676f2a66SJacob Faibussowitsch 
117676f2a66SJacob Faibussowitsch    when looking to insert "5" this routine will return *m = 6, but when looking to insert "6" it will ALSO return *m = 6.
118676f2a66SJacob Faibussowitsch   */
1199371c9d4SSatish Balay static inline PetscErrorCode PetscGallopSearchLeft_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m) {
120676f2a66SJacob Faibussowitsch   PetscInt last = l, k = 1, mid, cur = l + 1;
121676f2a66SJacob Faibussowitsch 
122676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
123676f2a66SJacob Faibussowitsch   *m = l;
1246bdcaf15SBarry Smith   PetscAssert(r >= l, PETSC_COMM_SELF, PETSC_ERR_PLIB, "r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchLeft", r, l);
1259371c9d4SSatish Balay   if ((*cmp)(x, arr + r * size, ctx) >= 0) {
1269371c9d4SSatish Balay     *m = r;
1279371c9d4SSatish Balay     PetscFunctionReturn(0);
1289371c9d4SSatish Balay   }
1294d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr) + l * size, ctx) < 0 || PetscUnlikely(!(r - l))) PetscFunctionReturn(0);
130676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
1319371c9d4SSatish Balay     if (PetscUnlikely(cur > r)) {
1329371c9d4SSatish Balay       cur = r;
1339371c9d4SSatish Balay       break;
1349371c9d4SSatish Balay     }
1354d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr) + cur * size, ctx) < 0) break;
136676f2a66SJacob Faibussowitsch     last = cur;
137676f2a66SJacob Faibussowitsch     cur += (k <<= 1) + 1;
138676f2a66SJacob Faibussowitsch     ++k;
139676f2a66SJacob Faibussowitsch   }
140676f2a66SJacob Faibussowitsch   /* standard binary search but take last 0 mid 0 cur 1 into account*/
141676f2a66SJacob Faibussowitsch   while (cur > last + 1) {
142676f2a66SJacob Faibussowitsch     mid = last + ((cur - last) >> 1);
1434d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr) + mid * size, ctx) < 0) {
144676f2a66SJacob Faibussowitsch       cur = mid;
145676f2a66SJacob Faibussowitsch     } else {
146676f2a66SJacob Faibussowitsch       last = mid;
147676f2a66SJacob Faibussowitsch     }
148676f2a66SJacob Faibussowitsch   }
149676f2a66SJacob Faibussowitsch   *m = cur;
150676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
151676f2a66SJacob Faibussowitsch }
152676f2a66SJacob Faibussowitsch 
153676f2a66SJacob 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]
154676f2a66SJacob Faibussowitsch  < x. Output also inclusive */
1559371c9d4SSatish Balay static inline PetscErrorCode PetscGallopSearchRight_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m) {
156676f2a66SJacob Faibussowitsch   PetscInt last = r, k = 1, mid, cur = r - 1;
157676f2a66SJacob Faibussowitsch 
158676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
159676f2a66SJacob Faibussowitsch   *m = r;
1606bdcaf15SBarry Smith   PetscAssert(r >= l, PETSC_COMM_SELF, PETSC_ERR_PLIB, "r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchRight", r, l);
1619371c9d4SSatish Balay   if ((*cmp)(x, arr + l * size, ctx) <= 0) {
1629371c9d4SSatish Balay     *m = l;
1639371c9d4SSatish Balay     PetscFunctionReturn(0);
1649371c9d4SSatish Balay   }
1654d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr) + r * size, ctx) > 0 || PetscUnlikely(!(r - l))) PetscFunctionReturn(0);
166676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
1679371c9d4SSatish Balay     if (PetscUnlikely(cur < l)) {
1689371c9d4SSatish Balay       cur = l;
1699371c9d4SSatish Balay       break;
1709371c9d4SSatish Balay     }
1714d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr) + cur * size, ctx) > 0) break;
172676f2a66SJacob Faibussowitsch     last = cur;
173676f2a66SJacob Faibussowitsch     cur -= (k <<= 1) + 1;
174676f2a66SJacob Faibussowitsch     ++k;
175676f2a66SJacob Faibussowitsch   }
176676f2a66SJacob Faibussowitsch   /* standard binary search but take last r-1 mid r-1 cur r-2 into account*/
177676f2a66SJacob Faibussowitsch   while (last > cur + 1) {
178676f2a66SJacob Faibussowitsch     mid = last - ((last - cur) >> 1);
1794d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr) + mid * size, ctx) > 0) {
180676f2a66SJacob Faibussowitsch       cur = mid;
181676f2a66SJacob Faibussowitsch     } else {
182676f2a66SJacob Faibussowitsch       last = mid;
183676f2a66SJacob Faibussowitsch     }
184676f2a66SJacob Faibussowitsch   }
185676f2a66SJacob Faibussowitsch   *m = cur;
186676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
187676f2a66SJacob Faibussowitsch }
188676f2a66SJacob Faibussowitsch 
189676f2a66SJacob Faibussowitsch /* Mergesort where size of left half <= size of right half, so mergesort is done left to right. Arr should be pointer to
190676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
191676f2a66SJacob Faibussowitsch  array */
1929371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortMergeLo_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) {
193676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
194676f2a66SJacob Faibussowitsch   const PetscInt llen = mid - left;
195676f2a66SJacob Faibussowitsch 
196676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
197676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, arr + (left * size), llen * size);
198676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
1994d3610e3SJacob Faibussowitsch     if ((*cmp)(tarr + (i * size), arr + (j * size), ctx) < 0) {
200676f2a66SJacob Faibussowitsch       Petsc_memcpy(arr + (k * size), tarr + (i * size), size);
201676f2a66SJacob Faibussowitsch       ++k;
202676f2a66SJacob Faibussowitsch       gallopright = 0;
203676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
204676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
205676f2a66SJacob Faibussowitsch         do {
206676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
207676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2089566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen - 1, arr + (j * size), &l1));
209676f2a66SJacob Faibussowitsch           diff1 = l1 - i;
210676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr + (k * size), tarr + (i * size), diff1 * size);
211676f2a66SJacob Faibussowitsch           k += diff1;
212676f2a66SJacob Faibussowitsch           i = l1;
213676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2149566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr + (i * size), &l2));
215676f2a66SJacob Faibussowitsch           diff2 = l2 - j;
216676f2a66SJacob Faibussowitsch           Petsc_memmove((arr) + k * size, (arr) + j * size, diff2 * size);
217676f2a66SJacob Faibussowitsch           k += diff2;
218676f2a66SJacob Faibussowitsch           j = l2;
219676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
220676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
221676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
222676f2a66SJacob Faibussowitsch       }
223676f2a66SJacob Faibussowitsch     } else {
224676f2a66SJacob Faibussowitsch       Petsc_memmove(arr + (k * size), arr + (j * size), size);
225676f2a66SJacob Faibussowitsch       ++k;
226676f2a66SJacob Faibussowitsch       gallopleft = 0;
227676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
228676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
229676f2a66SJacob Faibussowitsch         do {
230676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
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           /* search temp for right[j], can copy up to that of temp into arr immediately */
2389566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen - 1, arr + (j * size), &l1));
239676f2a66SJacob Faibussowitsch           diff1 = l1 - i;
240676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr + (k * size), tarr + (i * size), diff1 * size);
241676f2a66SJacob Faibussowitsch           k += diff1;
242676f2a66SJacob Faibussowitsch           i = l1;
243676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
244676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
245676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
246676f2a66SJacob Faibussowitsch       }
247676f2a66SJacob Faibussowitsch     }
248676f2a66SJacob Faibussowitsch   }
249ad540459SPierre Jolivet   if (i < llen) Petsc_memcpy(arr + (k * size), tarr + (i * size), (llen - i) * size);
250676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
251676f2a66SJacob Faibussowitsch }
252676f2a66SJacob Faibussowitsch 
2539371c9d4SSatish Balay 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) {
254676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
255676f2a66SJacob Faibussowitsch   const PetscInt llen = mid - left;
256676f2a66SJacob Faibussowitsch 
257676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
258676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, arr + (left * asize), llen * asize, btarr, barr + (left * bsize), llen * bsize);
259676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
2604d3610e3SJacob Faibussowitsch     if ((*cmp)(atarr + (i * asize), arr + (j * asize), ctx) < 0) {
261676f2a66SJacob Faibussowitsch       Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), asize, barr + (k * bsize), btarr + (i * bsize), bsize);
262676f2a66SJacob Faibussowitsch       ++k;
263676f2a66SJacob Faibussowitsch       gallopright = 0;
264676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
265676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
266676f2a66SJacob Faibussowitsch         do {
267676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
268676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2699566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen - 1, arr + (j * asize), &l1));
270676f2a66SJacob Faibussowitsch           diff1 = l1 - i;
271676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), diff1 * asize, barr + (k * bsize), btarr + (i * bsize), diff1 * bsize);
272676f2a66SJacob Faibussowitsch           k += diff1;
273676f2a66SJacob Faibussowitsch           i = l1;
274676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2759566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr + (i * asize), &l2));
276676f2a66SJacob Faibussowitsch           diff2 = l2 - j;
277676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr + (k * asize), arr + (j * asize), diff2 * asize, barr + (k * bsize), barr + (j * bsize), diff2 * bsize);
278676f2a66SJacob Faibussowitsch           k += diff2;
279676f2a66SJacob Faibussowitsch           j = l2;
280676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
281676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
282676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
283676f2a66SJacob Faibussowitsch       }
284676f2a66SJacob Faibussowitsch     } else {
285676f2a66SJacob Faibussowitsch       Petsc_memmove2(arr + (k * asize), arr + (j * asize), asize, barr + (k * bsize), barr + (j * bsize), bsize);
286676f2a66SJacob Faibussowitsch       ++k;
287676f2a66SJacob Faibussowitsch       gallopleft = 0;
288676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
289676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
290676f2a66SJacob Faibussowitsch         do {
291676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
292676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2939566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr + (i * asize), &l2));
294676f2a66SJacob Faibussowitsch           diff2 = l2 - j;
295676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr + (k * asize), arr + (j * asize), diff2 * asize, barr + (k * bsize), barr + (j * bsize), diff2 * bsize);
296676f2a66SJacob Faibussowitsch           k += diff2;
297676f2a66SJacob Faibussowitsch           j = l2;
298676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
2999566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen - 1, arr + (j * asize), &l1));
300676f2a66SJacob Faibussowitsch           diff1 = l1 - i;
301676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), diff1 * asize, barr + (k * bsize), btarr + (i * bsize), diff1 * bsize);
302676f2a66SJacob Faibussowitsch           k += diff1;
303676f2a66SJacob Faibussowitsch           i = l1;
304676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
305676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
306676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
307676f2a66SJacob Faibussowitsch       }
308676f2a66SJacob Faibussowitsch     }
309676f2a66SJacob Faibussowitsch   }
310ad540459SPierre Jolivet   if (i < llen) Petsc_memcpy2(arr + (k * asize), atarr + (i * asize), (llen - i) * asize, barr + (k * bsize), btarr + (i * bsize), (llen - i) * bsize);
311676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
312676f2a66SJacob Faibussowitsch }
313676f2a66SJacob Faibussowitsch 
314676f2a66SJacob Faibussowitsch /* Mergesort where size of right half < size of left half, so mergesort is done right to left. Arr should be pointer to
315676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
316676f2a66SJacob Faibussowitsch  array */
3179371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortMergeHi_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right) {
318676f2a66SJacob Faibussowitsch   PetscInt       i = right - mid, j = mid - 1, k = right, gallopleft = 0, gallopright = 0;
319676f2a66SJacob Faibussowitsch   const PetscInt rlen = right - mid + 1;
320676f2a66SJacob Faibussowitsch 
321676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
322676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, (arr) + mid * size, rlen * size);
323676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
3244d3610e3SJacob Faibussowitsch     if ((*cmp)((tarr) + i * size, (arr) + j * size, ctx) > 0) {
325676f2a66SJacob Faibussowitsch       Petsc_memcpy((arr) + k * size, (tarr) + i * size, size);
326676f2a66SJacob Faibussowitsch       --k;
327676f2a66SJacob Faibussowitsch       gallopleft = 0;
328676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
329676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
330676f2a66SJacob Faibussowitsch         do {
331676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
332676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3339566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr) + j * size, &l1));
334676f2a66SJacob Faibussowitsch           diff1 = i - l1;
335676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr) + (k - diff1 + 1) * size, (tarr) + (l1 + 1) * size, diff1 * size);
336676f2a66SJacob Faibussowitsch           k -= diff1;
337676f2a66SJacob Faibussowitsch           i = l1;
338676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3399566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr) + i * size, &l2));
340676f2a66SJacob Faibussowitsch           diff2 = j - l2;
341676f2a66SJacob Faibussowitsch           Petsc_memmove((arr) + (k - diff2 + 1) * size, (arr) + (l2 + 1) * size, diff2 * size);
342676f2a66SJacob Faibussowitsch           k -= diff2;
343676f2a66SJacob Faibussowitsch           j = l2;
344676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
345676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
346676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
347676f2a66SJacob Faibussowitsch       }
348676f2a66SJacob Faibussowitsch     } else {
349676f2a66SJacob Faibussowitsch       Petsc_memmove((arr) + k * size, (arr) + j * size, size);
350676f2a66SJacob Faibussowitsch       --k;
351676f2a66SJacob Faibussowitsch       gallopright = 0;
352676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
353676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
354676f2a66SJacob Faibussowitsch         do {
355676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
356676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3579566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr) + i * size, &l2));
358676f2a66SJacob Faibussowitsch           diff2 = j - l2;
359676f2a66SJacob Faibussowitsch           Petsc_memmove((arr) + (k - diff2 + 1) * size, (arr) + (l2 + 1) * size, diff2 * size);
360676f2a66SJacob Faibussowitsch           k -= diff2;
361676f2a66SJacob Faibussowitsch           j = l2;
362676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3639566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr) + j * size, &l1));
364676f2a66SJacob Faibussowitsch           diff1 = i - l1;
365676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr) + (k - diff1 + 1) * size, (tarr) + (l1 + 1) * size, diff1 * size);
366676f2a66SJacob Faibussowitsch           k -= diff1;
367676f2a66SJacob Faibussowitsch           i = l1;
368676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
369676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
370676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
371676f2a66SJacob Faibussowitsch       }
372676f2a66SJacob Faibussowitsch     }
373676f2a66SJacob Faibussowitsch   }
374ad540459SPierre Jolivet   if (i >= 0) Petsc_memcpy((arr) + left * size, tarr, (i + 1) * size);
375676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
376676f2a66SJacob Faibussowitsch }
377676f2a66SJacob Faibussowitsch 
3789371c9d4SSatish Balay 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) {
379676f2a66SJacob Faibussowitsch   PetscInt       i = right - mid, j = mid - 1, k = right, gallopleft = 0, gallopright = 0;
380676f2a66SJacob Faibussowitsch   const PetscInt rlen = right - mid + 1;
381676f2a66SJacob Faibussowitsch 
382676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
383676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, (arr) + mid * asize, rlen * asize, btarr, (barr) + mid * bsize, rlen * bsize);
384676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
3854d3610e3SJacob Faibussowitsch     if ((*cmp)((atarr) + i * asize, (arr) + j * asize, ctx) > 0) {
386676f2a66SJacob Faibussowitsch       Petsc_memcpy2((arr) + k * asize, (atarr) + i * asize, asize, (barr) + k * bsize, (btarr) + i * bsize, bsize);
387676f2a66SJacob Faibussowitsch       --k;
388676f2a66SJacob Faibussowitsch       gallopleft = 0;
389676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
390676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
391676f2a66SJacob Faibussowitsch         do {
392676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
393676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3949566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr) + j * asize, &l1));
395676f2a66SJacob Faibussowitsch           diff1 = i - l1;
396676f2a66SJacob 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);
397676f2a66SJacob Faibussowitsch           k -= diff1;
398676f2a66SJacob Faibussowitsch           i = l1;
399676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4009566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr) + i * asize, &l2));
401676f2a66SJacob Faibussowitsch           diff2 = j - l2;
402676f2a66SJacob 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);
403676f2a66SJacob Faibussowitsch           k -= diff2;
404676f2a66SJacob Faibussowitsch           j = l2;
405676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
406676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
407676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
408676f2a66SJacob Faibussowitsch       }
409676f2a66SJacob Faibussowitsch     } else {
410676f2a66SJacob Faibussowitsch       Petsc_memmove2((arr) + k * asize, (arr) + j * asize, asize, (barr) + k * bsize, (barr) + j * bsize, bsize);
411676f2a66SJacob Faibussowitsch       --k;
412676f2a66SJacob Faibussowitsch       gallopright = 0;
413676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
414676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
415676f2a66SJacob Faibussowitsch         do {
416676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
417676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4189566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr) + i * asize, &l2));
419676f2a66SJacob Faibussowitsch           diff2 = j - l2;
420676f2a66SJacob 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);
421676f2a66SJacob Faibussowitsch           k -= diff2;
422676f2a66SJacob Faibussowitsch           j = l2;
423676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
4249566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr) + j * asize, &l1));
425676f2a66SJacob Faibussowitsch           diff1 = i - l1;
426676f2a66SJacob 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);
427676f2a66SJacob Faibussowitsch           k -= diff1;
428676f2a66SJacob Faibussowitsch           i = l1;
429676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
430676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
431676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
432676f2a66SJacob Faibussowitsch       }
433676f2a66SJacob Faibussowitsch     }
434676f2a66SJacob Faibussowitsch   }
435ad540459SPierre Jolivet   if (i >= 0) Petsc_memcpy2((arr) + left * asize, atarr, (i + 1) * asize, (barr) + left * bsize, btarr, (i + 1) * bsize);
436676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
437676f2a66SJacob Faibussowitsch }
438676f2a66SJacob Faibussowitsch 
439676f2a66SJacob Faibussowitsch /* Left is inclusive lower bound of array slice, start is start location of unsorted section, right is inclusive upper
440676f2a66SJacob Faibussowitsch  bound of array slice. If unsure of where unsorted section starts or if entire length is unsorted pass start = left */
4419371c9d4SSatish Balay static inline PetscErrorCode PetscInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) {
442676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start + 1 : start;
443676f2a66SJacob Faibussowitsch 
444676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
445676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
446676f2a66SJacob Faibussowitsch     PetscInt j = i - 1;
447676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr + (i * size), size);
4484d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(tarr, (arr) + j * size, ctx) < 0)) {
449676f2a66SJacob Faibussowitsch       COPYSWAPPY(arr + (j + 1) * size, arr + j * size, tarr + size, size);
450676f2a66SJacob Faibussowitsch       --j;
451676f2a66SJacob Faibussowitsch     }
452676f2a66SJacob Faibussowitsch     Petsc_memcpy((arr) + (j + 1) * size, tarr, size);
453676f2a66SJacob Faibussowitsch   }
454676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
455676f2a66SJacob Faibussowitsch }
456676f2a66SJacob Faibussowitsch 
4579371c9d4SSatish Balay 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) {
458676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start + 1 : start;
459676f2a66SJacob Faibussowitsch 
460676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
461676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
462676f2a66SJacob Faibussowitsch     PetscInt j = i - 1;
463676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr + (i * asize), asize, btarr, barr + (i * bsize), bsize);
4644d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(atarr, arr + (j * asize), ctx) < 0)) {
465676f2a66SJacob Faibussowitsch       COPYSWAPPY2(arr + (j + 1) * asize, arr + j * asize, asize, barr + (j + 1) * bsize, barr + j * bsize, bsize, atarr + asize);
466676f2a66SJacob Faibussowitsch       --j;
467676f2a66SJacob Faibussowitsch     }
468676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr + (j + 1) * asize, atarr, asize, barr + (j + 1) * bsize, btarr, bsize);
469676f2a66SJacob Faibussowitsch   }
470676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
471676f2a66SJacob Faibussowitsch }
472676f2a66SJacob Faibussowitsch 
473676f2a66SJacob Faibussowitsch /* See PetscInsertionSort_Private */
4749371c9d4SSatish Balay static inline PetscErrorCode PetscBinaryInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right) {
475676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start + 1 : start;
476676f2a66SJacob Faibussowitsch 
477676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
478676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
479676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
480676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr + (i * size), size);
481676f2a66SJacob Faibussowitsch     do {
482676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
4834d3610e3SJacob Faibussowitsch       if ((*cmp)(tarr, arr + (m * size), ctx) < 0) {
484676f2a66SJacob Faibussowitsch         r = m;
485676f2a66SJacob Faibussowitsch       } else {
486676f2a66SJacob Faibussowitsch         l = m + 1;
487676f2a66SJacob Faibussowitsch       }
488676f2a66SJacob Faibussowitsch     } while (l < r);
489676f2a66SJacob Faibussowitsch     Petsc_memmove(arr + ((l + 1) * size), arr + (l * size), (i - l) * size);
490676f2a66SJacob Faibussowitsch     Petsc_memcpy(arr + (l * size), tarr, size);
491676f2a66SJacob Faibussowitsch   }
492676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
493676f2a66SJacob Faibussowitsch }
494676f2a66SJacob Faibussowitsch 
4959371c9d4SSatish Balay 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) {
496676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start + 1 : start;
497676f2a66SJacob Faibussowitsch 
498676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
499676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
500676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
501676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr + (i * asize), asize, btarr, barr + (i * bsize), bsize);
502676f2a66SJacob Faibussowitsch     do {
503676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
5044d3610e3SJacob Faibussowitsch       if ((*cmp)(atarr, arr + (m * asize), ctx) < 0) {
505676f2a66SJacob Faibussowitsch         r = m;
506676f2a66SJacob Faibussowitsch       } else {
507676f2a66SJacob Faibussowitsch         l = m + 1;
508676f2a66SJacob Faibussowitsch       }
509676f2a66SJacob Faibussowitsch     } while (l < r);
510676f2a66SJacob Faibussowitsch     Petsc_memmove2(arr + ((l + 1) * asize), arr + (l * asize), (i - l) * asize, barr + ((l + 1) * bsize), barr + (l * bsize), (i - l) * bsize);
511676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr + (l * asize), atarr, asize, barr + (l * bsize), btarr, bsize);
512676f2a66SJacob Faibussowitsch   }
513676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
514676f2a66SJacob Faibussowitsch }
515676f2a66SJacob Faibussowitsch 
516676f2a66SJacob Faibussowitsch typedef struct {
517676f2a66SJacob Faibussowitsch   PetscInt size;
518676f2a66SJacob Faibussowitsch   PetscInt start;
51989d292beSStefano Zampini #if defined(__CRAY_AARCH64) /* segfaults with Cray compilers for aarch64 on a64FX */
52089d292beSStefano Zampini } PetscTimSortStack;
52189d292beSStefano Zampini #else
522676f2a66SJacob Faibussowitsch } PetscTimSortStack PETSC_ATTRIBUTEALIGNED(2 * sizeof(PetscInt));
52389d292beSStefano Zampini #endif
524676f2a66SJacob Faibussowitsch 
525676f2a66SJacob Faibussowitsch typedef struct {
526676f2a66SJacob Faibussowitsch   char *ptr PETSC_ATTRIBUTEALIGNED(PETSC_MEMALIGN);
527676f2a66SJacob Faibussowitsch   size_t    size;
528676f2a66SJacob Faibussowitsch   size_t    maxsize;
529676f2a66SJacob Faibussowitsch } PetscTimSortBuffer;
530676f2a66SJacob Faibussowitsch 
5319371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortResizeBuffer_Private(PetscTimSortBuffer *buff, size_t newSize) {
532676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
533676f2a66SJacob Faibussowitsch   if (PetscLikely(newSize <= buff->size)) PetscFunctionReturn(0);
534676f2a66SJacob Faibussowitsch   {
535676f2a66SJacob Faibussowitsch     /* Can't be larger than n, there is merit to simply allocating buff to n to begin with */
536676f2a66SJacob Faibussowitsch     size_t newMax = PetscMin(newSize * newSize, buff->maxsize);
5379566063dSJacob Faibussowitsch     PetscCall(PetscFree(buff->ptr));
5389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(newMax, &buff->ptr));
539676f2a66SJacob Faibussowitsch     buff->size = newMax;
540676f2a66SJacob Faibussowitsch   }
541676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
542676f2a66SJacob Faibussowitsch }
543676f2a66SJacob Faibussowitsch 
5449371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortForceCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt stacksize) {
545676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
546676f2a66SJacob Faibussowitsch   for (; stacksize; --stacksize) {
547676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
5484d3610e3SJacob Faibussowitsch     if ((*cmp)(arr + (stack[stacksize].start - 1) * size, arr + (stack[stacksize].start) * size, ctx) > 0) {
549676f2a66SJacob Faibussowitsch       PetscInt l, m = stack[stacksize].start, r;
550676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
5519566063dSJacob Faibussowitsch       PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[stacksize - 1].start, stack[stacksize].start - 1, (arr) + (stack[stacksize].start) * size, &l));
552676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
5539566063dSJacob 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));
554676f2a66SJacob Faibussowitsch       if (m - l <= r - m) {
5559566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size));
5569566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
557676f2a66SJacob Faibussowitsch       } else {
5589566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size));
5599566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
560676f2a66SJacob Faibussowitsch       }
561676f2a66SJacob Faibussowitsch     }
562676f2a66SJacob Faibussowitsch     /* Update A with merge */
563676f2a66SJacob Faibussowitsch     stack[stacksize - 1].size += stack[stacksize].size;
564676f2a66SJacob Faibussowitsch   }
565676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
566676f2a66SJacob Faibussowitsch }
567676f2a66SJacob Faibussowitsch 
5689371c9d4SSatish Balay 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) {
569676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
570676f2a66SJacob Faibussowitsch   for (; stacksize; --stacksize) {
571676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
5724d3610e3SJacob Faibussowitsch     if ((*cmp)(arr + (stack[stacksize].start - 1) * asize, arr + (stack[stacksize].start) * asize, ctx) > 0) {
573676f2a66SJacob Faibussowitsch       PetscInt l, m = stack[stacksize].start, r;
574676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
5759566063dSJacob Faibussowitsch       PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[stacksize - 1].start, stack[stacksize].start - 1, (arr) + (stack[stacksize].start) * asize, &l));
576676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
5779566063dSJacob 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));
578676f2a66SJacob Faibussowitsch       if (m - l <= r - m) {
5799566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize));
5809566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize));
5819566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
582676f2a66SJacob Faibussowitsch       } else {
5839566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize));
5849566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize));
5859566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
586676f2a66SJacob Faibussowitsch       }
587676f2a66SJacob Faibussowitsch     }
588676f2a66SJacob Faibussowitsch     /* Update A with merge */
589676f2a66SJacob Faibussowitsch     stack[stacksize - 1].size += stack[stacksize].size;
590676f2a66SJacob Faibussowitsch   }
591676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
592676f2a66SJacob Faibussowitsch }
593676f2a66SJacob Faibussowitsch 
5949371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortMergeCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt *stacksize) {
595676f2a66SJacob Faibussowitsch   PetscInt i = *stacksize;
596676f2a66SJacob Faibussowitsch 
597676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
598676f2a66SJacob Faibussowitsch   while (i) {
599676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
600676f2a66SJacob Faibussowitsch 
601676f2a66SJacob Faibussowitsch     if (i == 1) {
602676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
603676f2a66SJacob Faibussowitsch       if (stack[i - 1].size < stack[i].size) {
604676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
6054d3610e3SJacob Faibussowitsch         if ((*cmp)(arr + (stack[i].start - 1) * size, arr + (stack[i].start) * size, ctx) > 0) {
606676f2a66SJacob Faibussowitsch           m = stack[i].start;
607676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
6089566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 1].start, stack[i].start - 1, (arr) + stack[i].start * size, &l));
609676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
6109566063dSJacob 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));
611676f2a66SJacob Faibussowitsch           if (m - l <= r - m) {
6129566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size));
6139566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
614676f2a66SJacob Faibussowitsch           } else {
6159566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size));
6169566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
617676f2a66SJacob Faibussowitsch           }
618676f2a66SJacob Faibussowitsch         }
619676f2a66SJacob Faibussowitsch         /* Update A with merge */
620676f2a66SJacob Faibussowitsch         stack[i - 1].size += stack[i].size;
621676f2a66SJacob Faibussowitsch         --i;
622676f2a66SJacob Faibussowitsch       }
623676f2a66SJacob Faibussowitsch     } else {
624676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
625676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
626676f2a66SJacob Faibussowitsch       if (stack[i - 2].size <= stack[i - 1].size + stack[i].size) {
627676f2a66SJacob Faibussowitsch         if (stack[i - 2].size < stack[i].size) {
628676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
6294d3610e3SJacob Faibussowitsch           if ((*cmp)(arr + (stack[i - 1].start - 1) * size, arr + (stack[i - 1].start) * size, ctx) > 0) {
630676f2a66SJacob Faibussowitsch             m = stack[i - 1].start;
631676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
6329566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 2].start, stack[i - 1].start - 1, (arr) + (stack[i - 1].start) * size, &l));
633676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
6349566063dSJacob 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));
635676f2a66SJacob Faibussowitsch             if (m - l <= r - m) {
6369566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size));
6379566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
638676f2a66SJacob Faibussowitsch             } else {
6399566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size));
6409566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
641676f2a66SJacob Faibussowitsch             }
642676f2a66SJacob Faibussowitsch           }
643676f2a66SJacob Faibussowitsch           /* Update A with merge */
644676f2a66SJacob Faibussowitsch           stack[i - 2].size += stack[i - 1].size;
645676f2a66SJacob Faibussowitsch           /* Push C up the stack */
646676f2a66SJacob Faibussowitsch           stack[i - 1].start = stack[i].start;
647676f2a66SJacob Faibussowitsch           stack[i - 1].size  = stack[i].size;
648676f2a66SJacob Faibussowitsch         } else {
649676f2a66SJacob Faibussowitsch         /* merge C into B */
650676f2a66SJacob Faibussowitsch         mergeBC:
651676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
6524d3610e3SJacob Faibussowitsch           if ((*cmp)(arr + (stack[i].start - 1) * size, arr + (stack[i].start) * size, ctx) > 0) {
653676f2a66SJacob Faibussowitsch             m = stack[i].start;
654676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
6559566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i - 1].start, stack[i].start - 1, arr + stack[i].start * size, &l));
656676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
6579566063dSJacob 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));
658676f2a66SJacob Faibussowitsch             if (m - l <= r - m) {
6599566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (m - l + 1) * size));
6609566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
661676f2a66SJacob Faibussowitsch             } else {
6629566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (r - m + 1) * size));
6639566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
664676f2a66SJacob Faibussowitsch             }
665676f2a66SJacob Faibussowitsch           }
666676f2a66SJacob Faibussowitsch           /* Update B with merge */
667676f2a66SJacob Faibussowitsch           stack[i - 1].size += stack[i].size;
668676f2a66SJacob Faibussowitsch         }
669676f2a66SJacob Faibussowitsch         --i;
670676f2a66SJacob Faibussowitsch       } else if (stack[i - 1].size <= stack[i].size) {
671676f2a66SJacob Faibussowitsch         /* merge C into B */
672676f2a66SJacob Faibussowitsch         goto mergeBC;
673676f2a66SJacob Faibussowitsch       }
674676f2a66SJacob Faibussowitsch     }
675676f2a66SJacob Faibussowitsch     if (itemp == i) break;
676676f2a66SJacob Faibussowitsch   }
677676f2a66SJacob Faibussowitsch   *stacksize = i;
678676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
679676f2a66SJacob Faibussowitsch }
680676f2a66SJacob Faibussowitsch 
6819371c9d4SSatish Balay 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) {
682676f2a66SJacob Faibussowitsch   PetscInt i = *stacksize;
683676f2a66SJacob Faibussowitsch 
684676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
685676f2a66SJacob Faibussowitsch   while (i) {
686676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
687676f2a66SJacob Faibussowitsch 
688676f2a66SJacob Faibussowitsch     if (i == 1) {
689676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
690676f2a66SJacob Faibussowitsch       if (stack[i - 1].size < stack[i].size) {
691676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
6924d3610e3SJacob Faibussowitsch         if ((*cmp)(arr + (stack[i].start - 1) * asize, arr + (stack[i].start) * asize, ctx) > 0) {
693676f2a66SJacob Faibussowitsch           m = stack[i].start;
694676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
6959566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 1].start, stack[i].start - 1, (arr) + stack[i].start * asize, &l));
696676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
6979566063dSJacob 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));
698676f2a66SJacob Faibussowitsch           if (m - l <= r - m) {
6999566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize));
7009566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize));
7019566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
702676f2a66SJacob Faibussowitsch           } else {
7039566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize));
7049566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize));
7059566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
706676f2a66SJacob Faibussowitsch           }
707676f2a66SJacob Faibussowitsch         }
708676f2a66SJacob Faibussowitsch         /* Update A with merge */
709676f2a66SJacob Faibussowitsch         stack[i - 1].size += stack[i].size;
710676f2a66SJacob Faibussowitsch         --i;
711676f2a66SJacob Faibussowitsch       }
712676f2a66SJacob Faibussowitsch     } else {
713676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
714676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
715676f2a66SJacob Faibussowitsch       if (stack[i - 2].size <= stack[i - 1].size + stack[i].size) {
716676f2a66SJacob Faibussowitsch         if (stack[i - 2].size < stack[i].size) {
717676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
7184d3610e3SJacob Faibussowitsch           if ((*cmp)(arr + (stack[i - 1].start - 1) * asize, arr + (stack[i - 1].start) * asize, ctx) > 0) {
719676f2a66SJacob Faibussowitsch             m = stack[i - 1].start;
720676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
7219566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 2].start, stack[i - 1].start - 1, (arr) + (stack[i - 1].start) * asize, &l));
722676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
7239566063dSJacob 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));
724676f2a66SJacob Faibussowitsch             if (m - l <= r - m) {
7259566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize));
7269566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize));
7279566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
728676f2a66SJacob Faibussowitsch             } else {
7299566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize));
7309566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize));
7319566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
732676f2a66SJacob Faibussowitsch             }
733676f2a66SJacob Faibussowitsch           }
734676f2a66SJacob Faibussowitsch           /* Update A with merge */
735676f2a66SJacob Faibussowitsch           stack[i - 2].size += stack[i - 1].size;
736676f2a66SJacob Faibussowitsch           /* Push C up the stack */
737676f2a66SJacob Faibussowitsch           stack[i - 1].start = stack[i].start;
738676f2a66SJacob Faibussowitsch           stack[i - 1].size  = stack[i].size;
739676f2a66SJacob Faibussowitsch         } else {
740676f2a66SJacob Faibussowitsch         /* merge C into B */
741676f2a66SJacob Faibussowitsch         mergeBC:
742676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
7434d3610e3SJacob Faibussowitsch           if ((*cmp)(arr + (stack[i].start - 1) * asize, arr + (stack[i].start) * asize, ctx) > 0) {
744676f2a66SJacob Faibussowitsch             m = stack[i].start;
745676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
7469566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i - 1].start, stack[i].start - 1, arr + stack[i].start * asize, &l));
747676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
7489566063dSJacob 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));
749676f2a66SJacob Faibussowitsch             if (m - l <= r - m) {
7509566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m - l + 1) * asize));
7519566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m - l + 1) * bsize));
7529566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
753676f2a66SJacob Faibussowitsch             } else {
7549566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r - m + 1) * asize));
7559566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r - m + 1) * bsize));
7569566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
757676f2a66SJacob Faibussowitsch             }
758676f2a66SJacob Faibussowitsch           }
759676f2a66SJacob Faibussowitsch           /* Update B with merge */
760676f2a66SJacob Faibussowitsch           stack[i - 1].size += stack[i].size;
761676f2a66SJacob Faibussowitsch         }
762676f2a66SJacob Faibussowitsch         --i;
763676f2a66SJacob Faibussowitsch       } else if (stack[i - 1].size <= stack[i].size) {
764676f2a66SJacob Faibussowitsch         /* merge C into B */
765676f2a66SJacob Faibussowitsch         goto mergeBC;
766676f2a66SJacob Faibussowitsch       }
767676f2a66SJacob Faibussowitsch     }
768676f2a66SJacob Faibussowitsch     if (itemp == i) break;
769676f2a66SJacob Faibussowitsch   }
770676f2a66SJacob Faibussowitsch   *stacksize = i;
771676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
772676f2a66SJacob Faibussowitsch }
773676f2a66SJacob Faibussowitsch 
774676f2a66SJacob Faibussowitsch /* March sequentially through the array building up a "run" of weakly increasing or strictly decreasing contiguous
775676f2a66SJacob Faibussowitsch  elements. Decreasing runs are reversed by swapping. If the run is less than minrun, artificially extend it via either
776676f2a66SJacob Faibussowitsch  binary insertion sort or regulat insertion sort */
7779371c9d4SSatish Balay static inline PetscErrorCode PetscTimSortBuildRun_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend) {
778676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart + minrun, n - 1);
779676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
780676f2a66SJacob Faibussowitsch 
781676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
7829371c9d4SSatish Balay   if (PetscUnlikely(runstart == n - 1)) {
7839371c9d4SSatish Balay     *runend = runstart;
7849371c9d4SSatish Balay     PetscFunctionReturn(0);
7859371c9d4SSatish Balay   }
786676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
7874d3610e3SJacob Faibussowitsch   if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) < 0) {
788676f2a66SJacob Faibussowitsch     ++ri;
789676f2a66SJacob Faibussowitsch     while (ri < n - 1) {
7904d3610e3SJacob Faibussowitsch       if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) >= 0) break;
791676f2a66SJacob Faibussowitsch       ++ri;
792676f2a66SJacob Faibussowitsch     }
793676f2a66SJacob Faibussowitsch     {
794676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
795ad540459SPierre Jolivet       for (; lo < hi; ++lo, --hi) COPYSWAPPY(arr + lo * size, arr + hi * size, tarr, size);
796676f2a66SJacob Faibussowitsch     }
797676f2a66SJacob Faibussowitsch   } else {
798676f2a66SJacob Faibussowitsch     ++ri;
799676f2a66SJacob Faibussowitsch     while (ri < n - 1) {
8004d3610e3SJacob Faibussowitsch       if ((*cmp)((arr) + (ri + 1) * size, (arr) + ri * size, ctx) < 0) break;
801676f2a66SJacob Faibussowitsch       ++ri;
802676f2a66SJacob Faibussowitsch     }
803676f2a66SJacob Faibussowitsch   }
804676f2a66SJacob Faibussowitsch   if (ri < re) {
805676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
806676f2a66SJacob Faibussowitsch      binary search */
807676f2a66SJacob Faibussowitsch     if (ri - runstart <= minrun >> 1) {
808676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8094d3610e3SJacob Faibussowitsch       PetscInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
810676f2a66SJacob Faibussowitsch     } else {
8114d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
812676f2a66SJacob Faibussowitsch     }
813676f2a66SJacob Faibussowitsch     *runend = re;
814676f2a66SJacob Faibussowitsch   } else *runend = ri;
815676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
816676f2a66SJacob Faibussowitsch }
817676f2a66SJacob Faibussowitsch 
8189371c9d4SSatish Balay 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) {
819676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart + minrun, n - 1);
820676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
821676f2a66SJacob Faibussowitsch 
822676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
8239371c9d4SSatish Balay   if (PetscUnlikely(runstart == n - 1)) {
8249371c9d4SSatish Balay     *runend = runstart;
8259371c9d4SSatish Balay     PetscFunctionReturn(0);
8269371c9d4SSatish Balay   }
827676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
8284d3610e3SJacob Faibussowitsch   if ((*cmp)((arr) + (ri + 1) * asize, arr + (ri * asize), ctx) < 0) {
829676f2a66SJacob Faibussowitsch     ++ri;
830676f2a66SJacob Faibussowitsch     while (ri < n - 1) {
8314d3610e3SJacob Faibussowitsch       if ((*cmp)((arr) + (ri + 1) * asize, (arr) + ri * asize, ctx) >= 0) break;
832676f2a66SJacob Faibussowitsch       ++ri;
833676f2a66SJacob Faibussowitsch     }
834676f2a66SJacob Faibussowitsch     {
835676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
836ad540459SPierre Jolivet       for (; lo < hi; ++lo, --hi) COPYSWAPPY2(arr + lo * asize, arr + hi * asize, asize, barr + lo * bsize, barr + hi * bsize, bsize, atarr);
837676f2a66SJacob Faibussowitsch     }
838676f2a66SJacob Faibussowitsch   } else {
839676f2a66SJacob Faibussowitsch     ++ri;
840676f2a66SJacob Faibussowitsch     while (ri < n - 1) {
8414d3610e3SJacob Faibussowitsch       if ((*cmp)((arr) + (ri + 1) * asize, (arr) + ri * asize, ctx) < 0) break;
842676f2a66SJacob Faibussowitsch       ++ri;
843676f2a66SJacob Faibussowitsch     }
844676f2a66SJacob Faibussowitsch   }
845676f2a66SJacob Faibussowitsch   if (ri < re) {
846676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
847676f2a66SJacob Faibussowitsch      binary search */
848676f2a66SJacob Faibussowitsch     if (ri - runstart <= minrun >> 1) {
849676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8504d3610e3SJacob Faibussowitsch       PetscInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
851676f2a66SJacob Faibussowitsch     } else {
8524d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
853676f2a66SJacob Faibussowitsch     }
854676f2a66SJacob Faibussowitsch     *runend = re;
855676f2a66SJacob Faibussowitsch   } else *runend = ri;
856676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
857676f2a66SJacob Faibussowitsch }
858676f2a66SJacob Faibussowitsch 
8594d3610e3SJacob Faibussowitsch /*@C
860676f2a66SJacob Faibussowitsch   PetscTimSort - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm.
861676f2a66SJacob Faibussowitsch 
862676f2a66SJacob Faibussowitsch   Not Collective
863676f2a66SJacob Faibussowitsch 
864676f2a66SJacob Faibussowitsch   Input Parameters:
865676f2a66SJacob Faibussowitsch + n    - number of values
866676f2a66SJacob Faibussowitsch . arr  - array to be sorted
867676f2a66SJacob Faibussowitsch . size - size in bytes of the datatype held in arr
8684d3610e3SJacob Faibussowitsch . cmp  - function pointer to comparison function
8694d3610e3SJacob Faibussowitsch - ctx  - optional context to be passed to comparison function, NULL if not needed
870676f2a66SJacob Faibussowitsch 
871676f2a66SJacob Faibussowitsch   Output Parameters:
872676f2a66SJacob Faibussowitsch . arr  - sorted array
873676f2a66SJacob Faibussowitsch 
8748da36234SJacob Faibussowitsch   Notes:
8758da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
8768da36234SJacob Faibussowitsch  sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims to
8778da36234SJacob Faibussowitsch  select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced arrays. To
8788da36234SJacob Faibussowitsch  do so it repeatedly triggers attempts throughout to merge adjacent runs.
8798da36234SJacob Faibussowitsch 
8808da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
8818da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
8828da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
8838da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
8848da36234SJacob Faibussowitsch   likely all contain similar data.
8858da36234SJacob Faibussowitsch 
886676f2a66SJacob Faibussowitsch   Sample usage:
887*811af0c4SBarry Smith   The comparison function must follow the `qsort()` comparison function paradigm, returning the sign of the difference
888676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
889676f2a66SJacob Faibussowitsch   may also
890676f2a66SJacob Faibussowitsch  change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if
891676f2a66SJacob Faibussowitsch  the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing
8928da36234SJacob Faibussowitsch   order
893676f2a66SJacob Faibussowitsch 
894676f2a66SJacob Faibussowitsch .vb
8954d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
896676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
897676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
898676f2a66SJacob Faibussowitsch   }
899676f2a66SJacob Faibussowitsch .ve
9004d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
9014d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
902676f2a66SJacob Faibussowitsch   Then pass the function
903676f2a66SJacob Faibussowitsch .vb
9044d3610e3SJacob Faibussowitsch   PetscTimSort(n, arr, sizeof(arr[0]), my_increasing_comparison_function, ctx)
905676f2a66SJacob Faibussowitsch .ve
906676f2a66SJacob Faibussowitsch 
9074d3610e3SJacob Faibussowitsch   Fortran Notes:
9084d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
9098da36234SJacob Faibussowitsch   returns result. For example
9104d3610e3SJacob Faibussowitsch .vb
9114d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
9124d3610e3SJacob Faibussowitsch    implicit none
9134d3610e3SJacob Faibussowitsch 
9144d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
9154d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
9164d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
9174d3610e3SJacob Faibussowitsch 
9184d3610e3SJacob Faibussowitsch    if (left < right) then
9194d3610e3SJacob Faibussowitsch      result = -1
9204d3610e3SJacob Faibussowitsch    else if (left == right) then
9214d3610e3SJacob Faibussowitsch      result = 0
9224d3610e3SJacob Faibussowitsch    else
9234d3610e3SJacob Faibussowitsch      result = 1
9244d3610e3SJacob Faibussowitsch    end if
9254d3610e3SJacob Faibussowitsch    return
9264d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
9274d3610e3SJacob Faibussowitsch .ve
9284d3610e3SJacob Faibussowitsch 
929676f2a66SJacob Faibussowitsch   References:
930606c0280SSatish Balay . * - Tim Peters. https://bugs.python.org/file4451/timsort.txt
931676f2a66SJacob Faibussowitsch 
932676f2a66SJacob Faibussowitsch   Level: developer
933676f2a66SJacob Faibussowitsch 
934db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscIntSortSemiOrdered()`, `PetscRealSortSemiOrdered()`, `PetscMPIIntSortSemiOrdered()`
9354d3610e3SJacob Faibussowitsch @*/
9369371c9d4SSatish Balay PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *, void *), void *ctx) {
937676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
938676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
939676f2a66SJacob Faibussowitsch   PetscTimSortBuffer buff;
940676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
941676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
942676f2a66SJacob Faibussowitsch 
943676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
944676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
945676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
946676f2a66SJacob Faibussowitsch   {
947676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
948676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
949676f2a66SJacob Faibussowitsch     while (t >= 64) {
950676f2a66SJacob Faibussowitsch       r |= t & 1;
951676f2a66SJacob Faibussowitsch       t >>= 1;
952676f2a66SJacob Faibussowitsch     }
953676f2a66SJacob Faibussowitsch     minrun = t + r;
954676f2a66SJacob Faibussowitsch   }
9554d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
9569566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun));
9574d3610e3SJacob Faibussowitsch     if (n < 64) {
9589566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "n %" PetscInt_FMT " < 64, consider using PetscSortInt() instead\n", n));
95908401ef6SPierre Jolivet     } else PetscCheck(!(minrun < 32) && !(minrun > 65), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Calculated minrun %" PetscInt_FMT " not in range (32,65)", minrun);
9604d3610e3SJacob Faibussowitsch   }
9619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t)minrun * size, &buff.ptr));
962676f2a66SJacob Faibussowitsch   buff.size         = (size_t)minrun * size;
963676f2a66SJacob Faibussowitsch   buff.maxsize      = (size_t)n * size;
964676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
965676f2a66SJacob Faibussowitsch   while (runstart < n) {
966676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
9679566063dSJacob Faibussowitsch     PetscCall(PetscTimSortBuildRun_Private((char *)arr, buff.ptr, size, cmp, ctx, n, minrun, runstart, &runend));
968676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
969676f2a66SJacob Faibussowitsch     runstack[stacksize].size  = runend - runstart + 1;
9709566063dSJacob Faibussowitsch     PetscCall(PetscTimSortMergeCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, &stacksize));
971676f2a66SJacob Faibussowitsch     ++stacksize;
972676f2a66SJacob Faibussowitsch     runstart = runend + 1;
973676f2a66SJacob Faibussowitsch   }
974676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
975676f2a66SJacob Faibussowitsch   --stacksize;
9769566063dSJacob Faibussowitsch   PetscCall(PetscTimSortForceCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, stacksize));
9779566063dSJacob Faibussowitsch   PetscCall(PetscFree(buff.ptr));
978676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
979676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
980676f2a66SJacob Faibussowitsch }
981676f2a66SJacob Faibussowitsch 
9824d3610e3SJacob Faibussowitsch /*@C
983676f2a66SJacob Faibussowitsch   PetscTimSortWithArray - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm and
984676f2a66SJacob Faibussowitsch   reorders a second array to match the first. The arrays need not be the same type.
985676f2a66SJacob Faibussowitsch 
986676f2a66SJacob Faibussowitsch   Not Collective
987676f2a66SJacob Faibussowitsch 
988676f2a66SJacob Faibussowitsch   Input Parameters:
989676f2a66SJacob Faibussowitsch + n     - number of values
990676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in arr
9916b867d5aSJose E. Roman . bsize - size in bytes of the datatype held in barr
9924d3610e3SJacob Faibussowitsch . cmp   - function pointer to comparison function
9934d3610e3SJacob Faibussowitsch - ctx   - optional context to be passed to comparison function, NULL if not needed
994676f2a66SJacob Faibussowitsch 
9956b867d5aSJose E. Roman   Input/Output Parameters:
9966b867d5aSJose E. Roman + arr  - array to be sorted, on output it is sorted
9976b867d5aSJose E. Roman - barr - array to be reordered, on output it is reordered
9988da36234SJacob Faibussowitsch 
9998da36234SJacob Faibussowitsch   Notes:
10008da36234SJacob 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
10018da36234SJacob Faibussowitsch   overlap.
10028da36234SJacob Faibussowitsch 
10038da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
10048da36234SJacob Faibussowitsch   sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims
10058da36234SJacob Faibussowitsch  to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced
10068da36234SJacob Faibussowitsch  arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs.
10078da36234SJacob Faibussowitsch 
10088da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
10098da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
10108da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
10118da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
10128da36234SJacob Faibussowitsch   likely all contain similar data.
1013676f2a66SJacob Faibussowitsch 
1014676f2a66SJacob Faibussowitsch   Sample usage:
1015*811af0c4SBarry Smith   The comparison function must follow the `qsort()` comparison function paradigm, returning the sign of the difference
1016676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
1017676f2a66SJacob Faibussowitsch   may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only
1018676f2a66SJacob Faibussowitsch   guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in
10198da36234SJacob Faibussowitsch   increasing order
1020676f2a66SJacob Faibussowitsch 
1021676f2a66SJacob Faibussowitsch .vb
10224d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
1023676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
1024676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
1025676f2a66SJacob Faibussowitsch   }
1026676f2a66SJacob Faibussowitsch .ve
10274d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
10284d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
1029676f2a66SJacob Faibussowitsch   Then pass the function
1030676f2a66SJacob Faibussowitsch .vb
10314d3610e3SJacob Faibussowitsch   PetscTimSortWithArray(n, arr, sizeof(arr[0]), barr, sizeof(barr[0]), my_increasing_comparison_function, ctx)
1032676f2a66SJacob Faibussowitsch .ve
1033676f2a66SJacob Faibussowitsch 
10344d3610e3SJacob Faibussowitsch   Fortran Notes:
10354d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
10368da36234SJacob Faibussowitsch   returns result. For example
10374d3610e3SJacob Faibussowitsch .vb
10384d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
10394d3610e3SJacob Faibussowitsch    implicit none
10404d3610e3SJacob Faibussowitsch 
10414d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
10424d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
10434d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
10444d3610e3SJacob Faibussowitsch 
10454d3610e3SJacob Faibussowitsch    if (left < right) then
10464d3610e3SJacob Faibussowitsch      result = -1
10474d3610e3SJacob Faibussowitsch    else if (left == right) then
10484d3610e3SJacob Faibussowitsch      result = 0
10494d3610e3SJacob Faibussowitsch    else
10504d3610e3SJacob Faibussowitsch      result = 1
10514d3610e3SJacob Faibussowitsch    end if
10524d3610e3SJacob Faibussowitsch    return
10534d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
10544d3610e3SJacob Faibussowitsch .ve
10554d3610e3SJacob Faibussowitsch 
1056676f2a66SJacob Faibussowitsch   References:
1057606c0280SSatish Balay . * - Tim Peters. https://bugs.python.org/file4451/timsort.txt
1058676f2a66SJacob Faibussowitsch 
1059676f2a66SJacob Faibussowitsch   Level: developer
1060676f2a66SJacob Faibussowitsch 
1061db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscIntSortSemiOrderedWithArray()`, `PetscRealSortSemiOrderedWithArrayInt()`, `PetscMPIIntSortSemiOrderedWithArray()`
10624d3610e3SJacob Faibussowitsch @*/
10639371c9d4SSatish Balay PetscErrorCode PetscTimSortWithArray(PetscInt n, void *arr, size_t asize, void *barr, size_t bsize, int (*cmp)(const void *, const void *, void *), void *ctx) {
1064676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
1065676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
1066676f2a66SJacob Faibussowitsch   PetscTimSortBuffer abuff, bbuff;
1067676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
1068676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
1069676f2a66SJacob Faibussowitsch 
1070676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1071676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
1072676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
1073676f2a66SJacob Faibussowitsch   {
1074676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
1075676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
1076676f2a66SJacob Faibussowitsch     while (t >= 64) {
1077676f2a66SJacob Faibussowitsch       r |= t & 1;
1078676f2a66SJacob Faibussowitsch       t >>= 1;
1079676f2a66SJacob Faibussowitsch     }
1080676f2a66SJacob Faibussowitsch     minrun = t + r;
1081676f2a66SJacob Faibussowitsch   }
10824d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
10839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun));
1084cc73adaaSBarry Smith     PetscCheck(n < 64 || (minrun >= 32 && minrun <= 65), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Calculated minrun %" PetscInt_FMT " not in range (32,65)", minrun);
10854d3610e3SJacob Faibussowitsch   }
10869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t)minrun * asize, &abuff.ptr));
1087676f2a66SJacob Faibussowitsch   abuff.size    = (size_t)minrun * asize;
1088676f2a66SJacob Faibussowitsch   abuff.maxsize = (size_t)n * asize;
10899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t)minrun * bsize, &bbuff.ptr));
1090676f2a66SJacob Faibussowitsch   bbuff.size        = (size_t)minrun * bsize;
1091676f2a66SJacob Faibussowitsch   bbuff.maxsize     = (size_t)n * bsize;
1092676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1093676f2a66SJacob Faibussowitsch   while (runstart < n) {
1094676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
10959566063dSJacob Faibussowitsch     PetscCall(PetscTimSortBuildRunWithArray_Private((char *)arr, abuff.ptr, asize, (char *)barr, bbuff.ptr, bsize, cmp, ctx, n, minrun, runstart, &runend));
1096676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
1097676f2a66SJacob Faibussowitsch     runstack[stacksize].size  = runend - runstart + 1;
10989566063dSJacob Faibussowitsch     PetscCall(PetscTimSortMergeCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, &stacksize));
1099676f2a66SJacob Faibussowitsch     ++stacksize;
1100676f2a66SJacob Faibussowitsch     runstart = runend + 1;
1101676f2a66SJacob Faibussowitsch   }
1102676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
1103676f2a66SJacob Faibussowitsch   --stacksize;
11049566063dSJacob Faibussowitsch   PetscCall(PetscTimSortForceCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, stacksize));
11059566063dSJacob Faibussowitsch   PetscCall(PetscFree(abuff.ptr));
11069566063dSJacob Faibussowitsch   PetscCall(PetscFree(bbuff.ptr));
1107676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1108676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1109676f2a66SJacob Faibussowitsch }
1110676f2a66SJacob Faibussowitsch 
1111676f2a66SJacob Faibussowitsch /*@
1112*811af0c4SBarry Smith    PetscIntSortSemiOrdered - Sorts an array of `PetscInt` in place in increasing order.
1113676f2a66SJacob Faibussowitsch 
1114676f2a66SJacob Faibussowitsch    Not Collective
1115676f2a66SJacob Faibussowitsch 
1116676f2a66SJacob Faibussowitsch    Input Parameters:
1117676f2a66SJacob Faibussowitsch +  n   - number of values
1118676f2a66SJacob Faibussowitsch -  arr - array of integers
1119676f2a66SJacob Faibussowitsch 
1120676f2a66SJacob Faibussowitsch    Output Parameters:
1121676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1122676f2a66SJacob Faibussowitsch 
1123676f2a66SJacob Faibussowitsch    Notes:
1124*811af0c4SBarry Smith    If the array is less than 64 entries long `PetscSortInt()` is automatically used.
1125676f2a66SJacob Faibussowitsch 
1126*811af0c4SBarry Smith    This function serves as an alternative to `PetscSortInt()`. While this function works for any array of integers it is
1127676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1128a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1129676f2a66SJacob Faibussowitsch 
1130676f2a66SJacob Faibussowitsch    Level: intermediate
1131676f2a66SJacob Faibussowitsch 
1132db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
1133676f2a66SJacob Faibussowitsch @*/
11349371c9d4SSatish Balay PetscErrorCode PetscIntSortSemiOrdered(PetscInt n, PetscInt arr[]) {
1135676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1136676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1137676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr, 2);
1138676f2a66SJacob Faibussowitsch   if (n < 64) {
11399566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(n, arr));
1140676f2a66SJacob Faibussowitsch   } else {
11419566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscInt), Compare_PetscInt_Private, NULL));
1142676f2a66SJacob Faibussowitsch   }
1143676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1144676f2a66SJacob Faibussowitsch }
1145676f2a66SJacob Faibussowitsch 
1146676f2a66SJacob Faibussowitsch /*@
1147*811af0c4SBarry Smith    PetscIntSortSemiOrderedWithArray - Sorts an array of `PetscInt` in place in increasing order and reorders a second
1148*811af0c4SBarry Smith    `PetscInt` array to match the first.
1149676f2a66SJacob Faibussowitsch 
1150676f2a66SJacob Faibussowitsch    Not Collective
1151676f2a66SJacob Faibussowitsch 
11526b867d5aSJose E. Roman    Input Parameter:
11536b867d5aSJose E. Roman .  n   - number of values
1154676f2a66SJacob Faibussowitsch 
11556b867d5aSJose E. Roman    Input/Output Parameters:
11566b867d5aSJose E. Roman +  arr1 - array of integers to be sorted, modified on output
11576b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1158676f2a66SJacob Faibussowitsch 
1159676f2a66SJacob Faibussowitsch    Notes:
1160676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1161676f2a66SJacob Faibussowitsch 
1162676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortIntWithArray(). While this function works for any array of integers it is
1163676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1164a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1165676f2a66SJacob Faibussowitsch 
1166676f2a66SJacob Faibussowitsch    Level: intermediate
1167676f2a66SJacob Faibussowitsch 
1168db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortIntWithArray()`, `PetscSortIntWithPermutation()`
1169676f2a66SJacob Faibussowitsch @*/
11709371c9d4SSatish Balay PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt n, PetscInt arr1[], PetscInt arr2[]) {
1171676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
117252e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1173676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1, 2);
1174676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2, 3);
11756e260675SJacob Faibussowitsch   /* cannot export out to PetscIntSortWithArray here since it isn't stable */
11769566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscInt), arr2, sizeof(PetscInt), Compare_PetscInt_Private, NULL));
1177676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1178676f2a66SJacob Faibussowitsch }
1179676f2a66SJacob Faibussowitsch 
1180676f2a66SJacob Faibussowitsch /*@
1181*811af0c4SBarry Smith    PetscMPIIntSortSemiOrdered - Sorts an array of `PetscMPIInt` in place in increasing order.
1182676f2a66SJacob Faibussowitsch 
1183676f2a66SJacob Faibussowitsch    Not Collective
1184676f2a66SJacob Faibussowitsch 
1185676f2a66SJacob Faibussowitsch    Input Parameters:
1186676f2a66SJacob Faibussowitsch +  n   - number of values
1187*811af0c4SBarry Smith -  arr - array of `PetscMPIInt`
1188676f2a66SJacob Faibussowitsch 
1189676f2a66SJacob Faibussowitsch    Output Parameters:
1190676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1191676f2a66SJacob Faibussowitsch 
1192676f2a66SJacob Faibussowitsch    Notes:
1193*811af0c4SBarry Smith    If the array is less than 64 entries long `PetscSortMPIInt()` is automatically used.
1194676f2a66SJacob Faibussowitsch 
1195*811af0c4SBarry Smith    This function serves as an alternative to `PetscSortMPIInt()`. While this function works for any array of `PetscMPIInt` 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 
1199676f2a66SJacob Faibussowitsch    Level: intermediate
1200676f2a66SJacob Faibussowitsch 
1201db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortMPIInt()`
1202676f2a66SJacob Faibussowitsch @*/
12039371c9d4SSatish Balay PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt n, PetscMPIInt arr[]) {
1204676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
120552e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1206676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr, 2);
1207676f2a66SJacob Faibussowitsch   if (n < 64) {
12089566063dSJacob Faibussowitsch     PetscCall(PetscSortMPIInt(n, arr));
1209676f2a66SJacob Faibussowitsch   } else {
12109566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL));
1211676f2a66SJacob Faibussowitsch   }
1212676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1213676f2a66SJacob Faibussowitsch }
1214676f2a66SJacob Faibussowitsch 
1215676f2a66SJacob Faibussowitsch /*@
1216*811af0c4SBarry Smith    PetscMPIIntSortSemiOrderedWithArray - Sorts an array of `PetscMPIInt` in place in increasing order and reorders a second `PetscMPIInt`
1217676f2a66SJacob Faibussowitsch    array to match the first.
1218676f2a66SJacob Faibussowitsch 
1219676f2a66SJacob Faibussowitsch    Not Collective
1220676f2a66SJacob Faibussowitsch 
12216b867d5aSJose E. Roman    Input Parameter:
12226b867d5aSJose E. Roman .  n   - number of values
1223676f2a66SJacob Faibussowitsch 
12246b867d5aSJose E. Roman    Input/Output Parameters:
12256b867d5aSJose E. Roman .  arr1 - array of integers to be sorted, modified on output
12266b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1227676f2a66SJacob Faibussowitsch 
1228676f2a66SJacob Faibussowitsch    Notes:
1229676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1230676f2a66SJacob Faibussowitsch 
1231*811af0c4SBarry Smith    This function serves as an alternative to `PetscSortMPIIntWithArray()`. While this function works for any array of integers it is
1232676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1233a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1234676f2a66SJacob Faibussowitsch 
1235676f2a66SJacob Faibussowitsch    Level: intermediate
1236676f2a66SJacob Faibussowitsch 
1237db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortMPIIntWithArray()`, `PetscSortMPIIntWithPermutation()`
1238676f2a66SJacob Faibussowitsch @*/
12399371c9d4SSatish Balay PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt n, PetscMPIInt arr1[], PetscMPIInt arr2[]) {
1240676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1241676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1242676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1, 2);
1243676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2, 3);
12446e260675SJacob Faibussowitsch   /* cannot export out to PetscMPIIntSortWithArray here since it isn't stable */
12459566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscMPIInt), arr2, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL));
1246676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1247676f2a66SJacob Faibussowitsch }
1248676f2a66SJacob Faibussowitsch 
1249676f2a66SJacob Faibussowitsch /*@
1250*811af0c4SBarry Smith    PetscRealSortSemiOrdered - Sorts an array of `PetscReal` in place in increasing order.
1251676f2a66SJacob Faibussowitsch 
1252676f2a66SJacob Faibussowitsch    Not Collective
1253676f2a66SJacob Faibussowitsch 
1254676f2a66SJacob Faibussowitsch    Input Parameters:
1255676f2a66SJacob Faibussowitsch +  n   - number of values
1256*811af0c4SBarry Smith -  arr - array of `PetscReal`
1257676f2a66SJacob Faibussowitsch 
1258676f2a66SJacob Faibussowitsch    Output Parameters:
1259676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1260676f2a66SJacob Faibussowitsch 
1261676f2a66SJacob Faibussowitsch    Notes:
1262*811af0c4SBarry Smith    If the array is less than 64 entries long `PetscSortReal()` is automatically used.
1263676f2a66SJacob Faibussowitsch 
1264*811af0c4SBarry Smith    This function serves as an alternative to `PetscSortReal()`. While this function works for any array of `PetscReal` it is
1265676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1266a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1267676f2a66SJacob Faibussowitsch 
1268676f2a66SJacob Faibussowitsch    Level: intermediate
1269676f2a66SJacob Faibussowitsch 
1270db781477SPatrick Sanan .seealso: `PetscTimSort()`, `PetscSortReal()`, `PetscSortRealWithPermutation()`
1271676f2a66SJacob Faibussowitsch @*/
12729371c9d4SSatish Balay PetscErrorCode PetscRealSortSemiOrdered(PetscInt n, PetscReal arr[]) {
1273676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1274676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1275676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr, 2);
1276676f2a66SJacob Faibussowitsch   if (n < 64) {
12779566063dSJacob Faibussowitsch     PetscCall(PetscSortReal(n, arr));
1278676f2a66SJacob Faibussowitsch   } else {
12799566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscReal), Compare_PetscReal_Private, NULL));
1280676f2a66SJacob Faibussowitsch   }
1281676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1282676f2a66SJacob Faibussowitsch }
1283676f2a66SJacob Faibussowitsch 
1284676f2a66SJacob Faibussowitsch /*@
1285*811af0c4SBarry Smith    PetscRealSortSemiOrderedWithArrayInt - Sorts an array of `PetscReal` in place in increasing order and reorders a second
1286*811af0c4SBarry Smith    array of `PetscInt` to match the first.
1287676f2a66SJacob Faibussowitsch 
1288676f2a66SJacob Faibussowitsch    Not Collective
1289676f2a66SJacob Faibussowitsch 
12906b867d5aSJose E. Roman    Input Parameter:
12916b867d5aSJose E. Roman .  n   - number of values
1292676f2a66SJacob Faibussowitsch 
12936b867d5aSJose E. Roman    Input/Output Parameters:
1294*811af0c4SBarry Smith .  arr1 - array of `PetscReal` to be sorted, modified on output
1295*811af0c4SBarry Smith -  arr2 - array of `PetscInt` to be reordered, modified on output
1296676f2a66SJacob Faibussowitsch 
1297676f2a66SJacob Faibussowitsch    Notes:
1298*811af0c4SBarry Smith    This function serves as an alternative to `PetscSortRealWithArray()`. While this function works for any array of `PetscReal` it is
1299676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1300a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1301676f2a66SJacob Faibussowitsch 
1302676f2a66SJacob Faibussowitsch    Level: intermediate
1303676f2a66SJacob Faibussowitsch 
1304db781477SPatrick Sanan .seealso: `PetscTimSortWithArray()`, `PetscSortRealWithArrayInt()`, `PetscSortRealWithPermutation()`
1305676f2a66SJacob Faibussowitsch @*/
13069371c9d4SSatish Balay PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt n, PetscReal arr1[], PetscInt arr2[]) {
1307676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1308676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1309676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr1, 2);
1310676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2, 3);
13116e260675SJacob Faibussowitsch   /* cannot export out to PetscRealSortWithArrayInt here since it isn't stable */
13129566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscReal), arr2, sizeof(PetscInt), Compare_PetscReal_Private, NULL));
1313676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1314676f2a66SJacob Faibussowitsch }
1315