xref: /petsc/src/sys/utils/sortso.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1676f2a66SJacob Faibussowitsch #include <petscsys.h>                /*I  "petscsys.h"  I*/
2676f2a66SJacob Faibussowitsch #include <petsc/private/petscimpl.h>
3676f2a66SJacob Faibussowitsch 
49fbee547SJacob Faibussowitsch static inline int Compare_PetscMPIInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx)
5676f2a66SJacob Faibussowitsch {
6676f2a66SJacob Faibussowitsch   PetscMPIInt l = *(PetscMPIInt *) left, r = *(PetscMPIInt *) right;
7676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
8676f2a66SJacob Faibussowitsch }
9676f2a66SJacob Faibussowitsch 
109fbee547SJacob Faibussowitsch static inline int Compare_PetscInt_Private(const void *left, const void *right, PETSC_UNUSED void *ctx)
11676f2a66SJacob Faibussowitsch {
12676f2a66SJacob Faibussowitsch   PetscInt l = *(PetscInt *) left, r = *(PetscInt *) right;
13676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
14676f2a66SJacob Faibussowitsch }
15676f2a66SJacob Faibussowitsch 
169fbee547SJacob Faibussowitsch static inline int Compare_PetscReal_Private(const void *left, const void *right, PETSC_UNUSED void *ctx)
17676f2a66SJacob Faibussowitsch {
18676f2a66SJacob Faibussowitsch   PetscReal l = *(PetscReal *) left, r = *(PetscReal *) right;
19676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
20676f2a66SJacob Faibussowitsch }
21676f2a66SJacob Faibussowitsch 
22676f2a66SJacob Faibussowitsch #define MIN_GALLOP_CONST_GLOBAL 8
23676f2a66SJacob Faibussowitsch static PetscInt MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
244d3610e3SJacob Faibussowitsch typedef int (*CompFunc)(const void *, const void *, void *);
25676f2a66SJacob Faibussowitsch 
26676f2a66SJacob Faibussowitsch /* Mostly to force clang uses the builtins, but can't hurt to have gcc compilers also do the same */
27676f2a66SJacob Faibussowitsch #if !defined (PETSC_USE_DEBUG) && defined(__GNUC__)
289fbee547SJacob Faibussowitsch static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size)
29676f2a66SJacob Faibussowitsch {
30676f2a66SJacob Faibussowitsch   __builtin_memcpy(t, b, size);
31676f2a66SJacob Faibussowitsch   __builtin_memmove(b, a, size);
32676f2a66SJacob Faibussowitsch   __builtin_memcpy(a, t, size);
33676f2a66SJacob Faibussowitsch   return;
34676f2a66SJacob Faibussowitsch }
35676f2a66SJacob Faibussowitsch 
369fbee547SJacob Faibussowitsch static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t)
37676f2a66SJacob Faibussowitsch {
386e260675SJacob Faibussowitsch   __builtin_memcpy(t, ar, asize);
396e260675SJacob Faibussowitsch   __builtin_memmove(ar, al, asize);
406e260675SJacob Faibussowitsch   __builtin_memcpy(al, t, asize);
416e260675SJacob Faibussowitsch   __builtin_memcpy(t, br, bsize);
426e260675SJacob Faibussowitsch   __builtin_memmove(br, bl, bsize);
436e260675SJacob Faibussowitsch   __builtin_memcpy(bl, t, bsize);
44676f2a66SJacob Faibussowitsch   return;
45676f2a66SJacob Faibussowitsch }
46676f2a66SJacob Faibussowitsch 
479fbee547SJacob Faibussowitsch static inline void Petsc_memcpy(char *dest, const char *src, size_t size)
48676f2a66SJacob Faibussowitsch {
49676f2a66SJacob Faibussowitsch   __builtin_memcpy(dest, src, size);
50676f2a66SJacob Faibussowitsch   return;
51676f2a66SJacob Faibussowitsch }
52676f2a66SJacob Faibussowitsch 
539fbee547SJacob Faibussowitsch static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
54676f2a66SJacob Faibussowitsch {
55676f2a66SJacob Faibussowitsch   __builtin_memcpy(adest, asrc, asize);
566e260675SJacob Faibussowitsch   __builtin_memcpy(bdest, bsrc, bsize);
57676f2a66SJacob Faibussowitsch   return;
58676f2a66SJacob Faibussowitsch }
59676f2a66SJacob Faibussowitsch 
609fbee547SJacob Faibussowitsch static inline void Petsc_memmove(char *dest, const char *src, size_t size)
61676f2a66SJacob Faibussowitsch {
62676f2a66SJacob Faibussowitsch   __builtin_memmove(dest, src, size);
63676f2a66SJacob Faibussowitsch   return;
64676f2a66SJacob Faibussowitsch }
65676f2a66SJacob Faibussowitsch 
669fbee547SJacob Faibussowitsch static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
67676f2a66SJacob Faibussowitsch {
68676f2a66SJacob Faibussowitsch   __builtin_memmove(adest, asrc, asize);
69676f2a66SJacob Faibussowitsch   __builtin_memmove(bdest, bsrc, bsize);
70676f2a66SJacob Faibussowitsch   return;
71676f2a66SJacob Faibussowitsch }
72676f2a66SJacob Faibussowitsch # else
739fbee547SJacob Faibussowitsch static inline void COPYSWAPPY(char *a, char *b, char *t, size_t size)
74676f2a66SJacob Faibussowitsch {
75676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(t, b, size));
779566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(b, a, size));
789566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(a, t, size));
79676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
80676f2a66SJacob Faibussowitsch }
81676f2a66SJacob Faibussowitsch 
829fbee547SJacob Faibussowitsch static inline void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t)
83676f2a66SJacob Faibussowitsch {
84676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(t, ar, asize));
869566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(ar, al, asize));
879566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(al, t, asize));
889566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(t, br, bsize));
899566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(br, bl, bsize));
909566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(bl, t, bsize));
91676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
92676f2a66SJacob Faibussowitsch }
93676f2a66SJacob Faibussowitsch 
949fbee547SJacob Faibussowitsch static inline void Petsc_memcpy(char *dest, const char *src, size_t size)
95676f2a66SJacob Faibussowitsch {
96676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(dest, src, size));
98676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
99676f2a66SJacob Faibussowitsch }
100676f2a66SJacob Faibussowitsch 
1019fbee547SJacob Faibussowitsch static inline void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
102676f2a66SJacob Faibussowitsch {
103676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(adest, asrc, asize));
1059566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemcpy(bdest, bsrc, bsize));
106676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
107676f2a66SJacob Faibussowitsch }
108676f2a66SJacob Faibussowitsch 
1099fbee547SJacob Faibussowitsch static inline void Petsc_memmove(char *dest, const char *src, size_t size)
110676f2a66SJacob Faibussowitsch {
111676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(dest, src, size));
113676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
114676f2a66SJacob Faibussowitsch }
115676f2a66SJacob Faibussowitsch 
1169fbee547SJacob Faibussowitsch static inline void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
117676f2a66SJacob Faibussowitsch {
118676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(adest, asrc, asize));
1209566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF,PetscMemmove(bdest, bsrc, bsize));
121676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
122676f2a66SJacob Faibussowitsch }
123676f2a66SJacob Faibussowitsch #endif
124676f2a66SJacob Faibussowitsch 
125676f2a66SJacob Faibussowitsch /* Start left look right. Looking for e.g. B[0] in A or mergelo. l inclusive, r inclusive. Returns first m such that arr[m] >
126676f2a66SJacob Faibussowitsch  x. Output also inclusive.
127676f2a66SJacob Faibussowitsch 
128676f2a66SJacob Faibussowitsch  NOTE: Both gallopsearch functions CANNOT distinguish between inserting AFTER the entire array vs inserting at entry n!! For example for an array:
129676f2a66SJacob Faibussowitsch 
130676f2a66SJacob Faibussowitsch    {0,1,2,3,4,5}
131676f2a66SJacob Faibussowitsch 
132676f2a66SJacob Faibussowitsch    when looking to insert "5" this routine will return *m = 6, but when looking to insert "6" it will ALSO return *m = 6.
133676f2a66SJacob Faibussowitsch   */
1349fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscGallopSearchLeft_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m)
135676f2a66SJacob Faibussowitsch {
136676f2a66SJacob Faibussowitsch   PetscInt last = l, k = 1, mid, cur = l+1;
137676f2a66SJacob Faibussowitsch 
138676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
139676f2a66SJacob Faibussowitsch   *m = l;
1406bdcaf15SBarry Smith   PetscAssert(r >= l,PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchLeft",r,l);
1414d3610e3SJacob Faibussowitsch   if ((*cmp)(x, arr+r*size, ctx) >= 0) {*m = r; PetscFunctionReturn(0);}
1424d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr)+l*size, ctx) < 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
143676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
144676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur > r)) {cur = r; break;}
1454d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size, ctx) < 0) break;
146676f2a66SJacob Faibussowitsch     last = cur;
147676f2a66SJacob Faibussowitsch     cur += (k <<= 1) + 1;
148676f2a66SJacob Faibussowitsch     ++k;
149676f2a66SJacob Faibussowitsch   }
150676f2a66SJacob Faibussowitsch   /* standard binary search but take last 0 mid 0 cur 1 into account*/
151676f2a66SJacob Faibussowitsch   while (cur > last + 1) {
152676f2a66SJacob Faibussowitsch     mid = last + ((cur - last) >> 1);
1534d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size, ctx) < 0) {
154676f2a66SJacob Faibussowitsch       cur = mid;
155676f2a66SJacob Faibussowitsch     } else {
156676f2a66SJacob Faibussowitsch       last = mid;
157676f2a66SJacob Faibussowitsch     }
158676f2a66SJacob Faibussowitsch   }
159676f2a66SJacob Faibussowitsch   *m = cur;
160676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
161676f2a66SJacob Faibussowitsch }
162676f2a66SJacob Faibussowitsch 
163676f2a66SJacob 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]
164676f2a66SJacob Faibussowitsch  < x. Output also inclusive */
1659fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscGallopSearchRight_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m)
166676f2a66SJacob Faibussowitsch {
167676f2a66SJacob Faibussowitsch   PetscInt last = r, k = 1, mid, cur = r-1;
168676f2a66SJacob Faibussowitsch 
169676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
170676f2a66SJacob Faibussowitsch   *m = r;
1716bdcaf15SBarry Smith   PetscAssert(r >= l,PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchRight",r,l);
1724d3610e3SJacob Faibussowitsch   if ((*cmp)(x, arr+l*size, ctx) <= 0) {*m = l; PetscFunctionReturn(0);}
1734d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr)+r*size, ctx) > 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
174676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
175676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur < l)) {cur = l; break;}
1764d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size, ctx) > 0) break;
177676f2a66SJacob Faibussowitsch     last = cur;
178676f2a66SJacob Faibussowitsch     cur -= (k <<= 1) + 1;
179676f2a66SJacob Faibussowitsch     ++k;
180676f2a66SJacob Faibussowitsch   }
181676f2a66SJacob Faibussowitsch   /* standard binary search but take last r-1 mid r-1 cur r-2 into account*/
182676f2a66SJacob Faibussowitsch   while (last > cur + 1) {
183676f2a66SJacob Faibussowitsch     mid = last - ((last - cur) >> 1);
1844d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size, ctx) > 0) {
185676f2a66SJacob Faibussowitsch       cur = mid;
186676f2a66SJacob Faibussowitsch     } else {
187676f2a66SJacob Faibussowitsch       last = mid;
188676f2a66SJacob Faibussowitsch     }
189676f2a66SJacob Faibussowitsch   }
190676f2a66SJacob Faibussowitsch   *m = cur;
191676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
192676f2a66SJacob Faibussowitsch }
193676f2a66SJacob Faibussowitsch 
194676f2a66SJacob Faibussowitsch /* Mergesort where size of left half <= size of right half, so mergesort is done left to right. Arr should be pointer to
195676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
196676f2a66SJacob Faibussowitsch  array */
1979fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeLo_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
198676f2a66SJacob Faibussowitsch {
199676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
200676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
201676f2a66SJacob Faibussowitsch 
202676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
203676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, arr+(left*size), llen*size);
204676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
2054d3610e3SJacob Faibussowitsch     if ((*cmp)(tarr+(i*size), arr+(j*size), ctx) < 0) {
206676f2a66SJacob Faibussowitsch       Petsc_memcpy(arr+(k*size), tarr+(i*size), size);
207676f2a66SJacob Faibussowitsch       ++k;
208676f2a66SJacob Faibussowitsch       gallopright = 0;
209676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
210676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
211676f2a66SJacob Faibussowitsch         do {
212676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
213676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2149566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen-1, arr+(j*size), &l1));
215676f2a66SJacob Faibussowitsch           diff1 = l1-i;
216676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
217676f2a66SJacob Faibussowitsch           k += diff1;
218676f2a66SJacob Faibussowitsch           i = l1;
219676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2209566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr+(i*size), &l2));
221676f2a66SJacob Faibussowitsch           diff2 = l2-j;
222676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+k*size, (arr)+j*size, diff2*size);
223676f2a66SJacob Faibussowitsch           k += diff2;
224676f2a66SJacob Faibussowitsch           j = l2;
225676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
226676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
227676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
228676f2a66SJacob Faibussowitsch       }
229676f2a66SJacob Faibussowitsch     } else {
230676f2a66SJacob Faibussowitsch       Petsc_memmove(arr+(k*size), arr+(j*size), size);
231676f2a66SJacob Faibussowitsch       ++k;
232676f2a66SJacob Faibussowitsch       gallopleft = 0;
233676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
234676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
235676f2a66SJacob Faibussowitsch         do {
236676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
237676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2389566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr+(i*size), &l2));
239676f2a66SJacob Faibussowitsch           diff2 = l2-j;
240676f2a66SJacob Faibussowitsch           Petsc_memmove(arr+(k*size), arr+(j*size), diff2*size);
241676f2a66SJacob Faibussowitsch           k += diff2;
242676f2a66SJacob Faibussowitsch           j = l2;
243676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
2449566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen-1, arr+(j*size), &l1));
245676f2a66SJacob Faibussowitsch           diff1 = l1-i;
246676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
247676f2a66SJacob Faibussowitsch           k += diff1;
248676f2a66SJacob Faibussowitsch           i = l1;
249676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
250676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
251676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
252676f2a66SJacob Faibussowitsch       }
253676f2a66SJacob Faibussowitsch     }
254676f2a66SJacob Faibussowitsch   }
255676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy(arr+(k*size), tarr+(i*size), (llen-i)*size);}
256676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
257676f2a66SJacob Faibussowitsch }
258676f2a66SJacob Faibussowitsch 
2599fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeLoWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
260676f2a66SJacob Faibussowitsch {
261676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
262676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
263676f2a66SJacob Faibussowitsch 
264676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
265676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, arr+(left*asize), llen*asize, btarr, barr+(left*bsize), llen*bsize);
266676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
2674d3610e3SJacob Faibussowitsch     if ((*cmp)(atarr+(i*asize), arr+(j*asize), ctx) < 0) {
268676f2a66SJacob Faibussowitsch       Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), asize, barr+(k*bsize), btarr+(i*bsize), bsize);
269676f2a66SJacob Faibussowitsch       ++k;
270676f2a66SJacob Faibussowitsch       gallopright = 0;
271676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
272676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
273676f2a66SJacob Faibussowitsch         do {
274676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
275676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2769566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen-1, arr+(j*asize), &l1));
277676f2a66SJacob Faibussowitsch           diff1 = l1-i;
278676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
279676f2a66SJacob Faibussowitsch           k += diff1;
280676f2a66SJacob Faibussowitsch           i = l1;
281676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2829566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr+(i*asize), &l2));
283676f2a66SJacob Faibussowitsch           diff2 = l2-j;
284676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
285676f2a66SJacob Faibussowitsch           k += diff2;
286676f2a66SJacob Faibussowitsch           j = l2;
287676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
288676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
289676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
290676f2a66SJacob Faibussowitsch       }
291676f2a66SJacob Faibussowitsch     } else {
292676f2a66SJacob Faibussowitsch       Petsc_memmove2(arr+(k*asize), arr+(j*asize), asize, barr+(k*bsize), barr+(j*bsize), bsize);
293676f2a66SJacob Faibussowitsch       ++k;
294676f2a66SJacob Faibussowitsch       gallopleft = 0;
295676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
296676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
297676f2a66SJacob Faibussowitsch         do {
298676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
299676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
3009566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr+(i*asize), &l2));
301676f2a66SJacob Faibussowitsch           diff2 = l2-j;
302676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
303676f2a66SJacob Faibussowitsch           k += diff2;
304676f2a66SJacob Faibussowitsch           j = l2;
305676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
3069566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen-1, arr+(j*asize), &l1));
307676f2a66SJacob Faibussowitsch           diff1 = l1-i;
308676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
309676f2a66SJacob Faibussowitsch           k += diff1;
310676f2a66SJacob Faibussowitsch           i = l1;
311676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
312676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
313676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
314676f2a66SJacob Faibussowitsch       }
315676f2a66SJacob Faibussowitsch     }
316676f2a66SJacob Faibussowitsch   }
317676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), (llen-i)*asize, barr+(k*bsize), btarr+(i*bsize), (llen-i)*bsize);}
318676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
319676f2a66SJacob Faibussowitsch }
320676f2a66SJacob Faibussowitsch 
321676f2a66SJacob Faibussowitsch /* Mergesort where size of right half < size of left half, so mergesort is done right to left. Arr should be pointer to
322676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
323676f2a66SJacob Faibussowitsch  array */
3249fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeHi_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
325676f2a66SJacob Faibussowitsch {
326676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
327676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
328676f2a66SJacob Faibussowitsch 
329676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
330676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, (arr)+mid*size, rlen*size);
331676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
3324d3610e3SJacob Faibussowitsch     if ((*cmp)((tarr)+i*size, (arr)+j*size, ctx) > 0) {
333676f2a66SJacob Faibussowitsch       Petsc_memcpy((arr)+k*size, (tarr)+i*size, size);
334676f2a66SJacob Faibussowitsch       --k;
335676f2a66SJacob Faibussowitsch       gallopleft = 0;
336676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
337676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
338676f2a66SJacob Faibussowitsch         do {
339676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
340676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3419566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr)+j*size, &l1));
342676f2a66SJacob Faibussowitsch           diff1 = i-l1;
343676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
344676f2a66SJacob Faibussowitsch           k -= diff1;
345676f2a66SJacob Faibussowitsch           i = l1;
346676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3479566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr)+i*size, &l2));
348676f2a66SJacob Faibussowitsch           diff2 = j-l2;
349676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
350676f2a66SJacob Faibussowitsch           k -= diff2;
351676f2a66SJacob Faibussowitsch           j = l2;
352676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
353676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
354676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
355676f2a66SJacob Faibussowitsch       }
356676f2a66SJacob Faibussowitsch     } else {
357676f2a66SJacob Faibussowitsch       Petsc_memmove((arr)+k*size, (arr)+j*size, size);
358676f2a66SJacob Faibussowitsch       --k;
359676f2a66SJacob Faibussowitsch       gallopright = 0;
360676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
361676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
362676f2a66SJacob Faibussowitsch         do {
363676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
364676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3659566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr)+i*size, &l2));
366676f2a66SJacob Faibussowitsch           diff2 = j-l2;
367676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
368676f2a66SJacob Faibussowitsch           k -= diff2;
369676f2a66SJacob Faibussowitsch           j = l2;
370676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3719566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr)+j*size, &l1));
372676f2a66SJacob Faibussowitsch           diff1 = i-l1;
373676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
374676f2a66SJacob Faibussowitsch           k -= diff1;
375676f2a66SJacob Faibussowitsch           i = l1;
376676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
377676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
378676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
379676f2a66SJacob Faibussowitsch       }
380676f2a66SJacob Faibussowitsch     }
381676f2a66SJacob Faibussowitsch   }
382676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy((arr)+left*size, tarr, (i+1)*size);}
383676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
384676f2a66SJacob Faibussowitsch }
385676f2a66SJacob Faibussowitsch 
3869fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeHiWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
387676f2a66SJacob Faibussowitsch {
388676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
389676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
390676f2a66SJacob Faibussowitsch 
391676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
392676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, (arr)+mid*asize, rlen*asize, btarr, (barr)+mid*bsize, rlen*bsize);
393676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
3944d3610e3SJacob Faibussowitsch     if ((*cmp)((atarr)+i*asize, (arr)+j*asize, ctx) > 0) {
395676f2a66SJacob Faibussowitsch       Petsc_memcpy2((arr)+k*asize, (atarr)+i*asize, asize, (barr)+k*bsize, (btarr)+i*bsize, bsize);
396676f2a66SJacob Faibussowitsch       --k;
397676f2a66SJacob Faibussowitsch       gallopleft = 0;
398676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
399676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
400676f2a66SJacob Faibussowitsch         do {
401676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
402676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
4039566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr)+j*asize, &l1));
404676f2a66SJacob Faibussowitsch           diff1 = i-l1;
405676f2a66SJacob 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);
406676f2a66SJacob Faibussowitsch           k -= diff1;
407676f2a66SJacob Faibussowitsch           i = l1;
408676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4099566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr)+i*asize, &l2));
410676f2a66SJacob Faibussowitsch           diff2 = j-l2;
411676f2a66SJacob 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);
412676f2a66SJacob Faibussowitsch           k -= diff2;
413676f2a66SJacob Faibussowitsch           j = l2;
414676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
415676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
416676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
417676f2a66SJacob Faibussowitsch       }
418676f2a66SJacob Faibussowitsch     } else {
419676f2a66SJacob Faibussowitsch       Petsc_memmove2((arr)+k*asize, (arr)+j*asize, asize, (barr)+k*bsize, (barr)+j*bsize, bsize);
420676f2a66SJacob Faibussowitsch       --k;
421676f2a66SJacob Faibussowitsch       gallopright = 0;
422676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
423676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
424676f2a66SJacob Faibussowitsch         do {
425676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
426676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4279566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr)+i*asize, &l2));
428676f2a66SJacob Faibussowitsch           diff2 = j-l2;
429676f2a66SJacob 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);
430676f2a66SJacob Faibussowitsch           k -= diff2;
431676f2a66SJacob Faibussowitsch           j = l2;
432676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
4339566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr)+j*asize, &l1));
434676f2a66SJacob Faibussowitsch           diff1 = i-l1;
435676f2a66SJacob 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);
436676f2a66SJacob Faibussowitsch           k -= diff1;
437676f2a66SJacob Faibussowitsch           i = l1;
438676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
439676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
440676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
441676f2a66SJacob Faibussowitsch       }
442676f2a66SJacob Faibussowitsch     }
443676f2a66SJacob Faibussowitsch   }
444676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy2((arr)+left*asize, atarr, (i+1)*asize, (barr)+left*bsize, btarr, (i+1)*bsize);}
445676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
446676f2a66SJacob Faibussowitsch }
447676f2a66SJacob Faibussowitsch 
448676f2a66SJacob Faibussowitsch /* Left is inclusive lower bound of array slice, start is start location of unsorted section, right is inclusive upper
449676f2a66SJacob Faibussowitsch  bound of array slice. If unsure of where unsorted section starts or if entire length is unsorted pass start = left */
4509fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
451676f2a66SJacob Faibussowitsch {
452676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
453676f2a66SJacob Faibussowitsch 
454676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
455676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
456676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
457676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
4584d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(tarr, (arr)+j*size, ctx) < 0)) {
459676f2a66SJacob Faibussowitsch       COPYSWAPPY(arr+(j+1)*size, arr+j*size, tarr+size, size);
460676f2a66SJacob Faibussowitsch       --j;
461676f2a66SJacob Faibussowitsch     }
462676f2a66SJacob Faibussowitsch     Petsc_memcpy((arr)+(j+1)*size, tarr, size);
463676f2a66SJacob Faibussowitsch   }
464676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
465676f2a66SJacob Faibussowitsch }
466676f2a66SJacob Faibussowitsch 
4679fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
468676f2a66SJacob Faibussowitsch {
469676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
470676f2a66SJacob Faibussowitsch 
471676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
472676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
473676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
474676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
4754d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(atarr, arr+(j*asize), ctx) < 0)) {
476676f2a66SJacob Faibussowitsch       COPYSWAPPY2(arr+(j+1)*asize, arr+j*asize, asize, barr+(j+1)*bsize, barr+j*bsize, bsize, atarr+asize);
477676f2a66SJacob Faibussowitsch       --j;
478676f2a66SJacob Faibussowitsch     }
479676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(j+1)*asize, atarr, asize, barr+(j+1)*bsize, btarr, bsize);
480676f2a66SJacob Faibussowitsch   }
481676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
482676f2a66SJacob Faibussowitsch }
483676f2a66SJacob Faibussowitsch 
484676f2a66SJacob Faibussowitsch /* See PetscInsertionSort_Private */
4859fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscBinaryInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
486676f2a66SJacob Faibussowitsch {
487676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
488676f2a66SJacob Faibussowitsch 
489676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
490676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
491676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
492676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
493676f2a66SJacob Faibussowitsch     do {
494676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
4954d3610e3SJacob Faibussowitsch       if ((*cmp)(tarr, arr+(m*size), ctx) < 0) {
496676f2a66SJacob Faibussowitsch         r = m;
497676f2a66SJacob Faibussowitsch       } else {
498676f2a66SJacob Faibussowitsch         l = m + 1;
499676f2a66SJacob Faibussowitsch       }
500676f2a66SJacob Faibussowitsch     } while (l < r);
501676f2a66SJacob Faibussowitsch     Petsc_memmove(arr+((l+1)*size), arr+(l*size), (i-l)*size);
502676f2a66SJacob Faibussowitsch     Petsc_memcpy(arr+(l*size), tarr, size);
503676f2a66SJacob Faibussowitsch   }
504676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
505676f2a66SJacob Faibussowitsch }
506676f2a66SJacob Faibussowitsch 
5079fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscBinaryInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
508676f2a66SJacob Faibussowitsch {
509676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
510676f2a66SJacob Faibussowitsch 
511676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
512676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
513676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
514676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
515676f2a66SJacob Faibussowitsch     do {
516676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
5174d3610e3SJacob Faibussowitsch       if ((*cmp)(atarr, arr+(m*asize), ctx) < 0) {
518676f2a66SJacob Faibussowitsch         r = m;
519676f2a66SJacob Faibussowitsch       } else {
520676f2a66SJacob Faibussowitsch         l = m + 1;
521676f2a66SJacob Faibussowitsch       }
522676f2a66SJacob Faibussowitsch     } while (l < r);
523676f2a66SJacob Faibussowitsch     Petsc_memmove2(arr+((l+1)*asize), arr+(l*asize), (i-l)*asize, barr+((l+1)*bsize), barr+(l*bsize), (i-l)*bsize);
524676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(l*asize), atarr, asize, barr+(l*bsize), btarr, bsize);
525676f2a66SJacob Faibussowitsch   }
526676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
527676f2a66SJacob Faibussowitsch }
528676f2a66SJacob Faibussowitsch 
529676f2a66SJacob Faibussowitsch typedef struct {
530676f2a66SJacob Faibussowitsch   PetscInt size;
531676f2a66SJacob Faibussowitsch   PetscInt start;
53289d292beSStefano Zampini #if defined(__CRAY_AARCH64) /* segfaults with Cray compilers for aarch64 on a64FX */
53389d292beSStefano Zampini } PetscTimSortStack;
53489d292beSStefano Zampini #else
535676f2a66SJacob Faibussowitsch } PetscTimSortStack PETSC_ATTRIBUTEALIGNED(2*sizeof(PetscInt));
53689d292beSStefano Zampini #endif
537676f2a66SJacob Faibussowitsch 
538676f2a66SJacob Faibussowitsch typedef struct {
539676f2a66SJacob Faibussowitsch   char   *ptr PETSC_ATTRIBUTEALIGNED(PETSC_MEMALIGN);
540676f2a66SJacob Faibussowitsch   size_t size;
541676f2a66SJacob Faibussowitsch   size_t maxsize;
542676f2a66SJacob Faibussowitsch } PetscTimSortBuffer;
543676f2a66SJacob Faibussowitsch 
5449fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortResizeBuffer_Private(PetscTimSortBuffer *buff, size_t newSize)
545676f2a66SJacob Faibussowitsch {
546676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
547676f2a66SJacob Faibussowitsch   if (PetscLikely(newSize <= buff->size)) PetscFunctionReturn(0);
548676f2a66SJacob Faibussowitsch   {
549676f2a66SJacob Faibussowitsch     /* Can't be larger than n, there is merit to simply allocating buff to n to begin with */
550676f2a66SJacob Faibussowitsch     size_t         newMax = PetscMin(newSize*newSize, buff->maxsize);
5519566063dSJacob Faibussowitsch     PetscCall(PetscFree(buff->ptr));
5529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(newMax, &buff->ptr));
553676f2a66SJacob Faibussowitsch     buff->size = newMax;
554676f2a66SJacob Faibussowitsch   }
555676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
556676f2a66SJacob Faibussowitsch }
557676f2a66SJacob Faibussowitsch 
5589fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortForceCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt stacksize)
559676f2a66SJacob Faibussowitsch {
560676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
561676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
562676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
5634d3610e3SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*size, arr+(stack[stacksize].start)*size, ctx) > 0) {
564676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
565676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
5669566063dSJacob Faibussowitsch       PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*size, &l));
567676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
5689566063dSJacob 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));
569676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
5709566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size));
5719566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
572676f2a66SJacob Faibussowitsch       } else {
5739566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size));
5749566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
575676f2a66SJacob Faibussowitsch       }
576676f2a66SJacob Faibussowitsch     }
577676f2a66SJacob Faibussowitsch     /* Update A with merge */
578676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
579676f2a66SJacob Faibussowitsch   }
580676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
581676f2a66SJacob Faibussowitsch }
582676f2a66SJacob Faibussowitsch 
5839fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortForceCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, void *ctx, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt stacksize)
584676f2a66SJacob Faibussowitsch {
585676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
586676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
587676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
5884d3610e3SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*asize, arr+(stack[stacksize].start)*asize, ctx) > 0) {
589676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
590676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
5919566063dSJacob Faibussowitsch       PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*asize, &l));
592676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
5939566063dSJacob 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));
594676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
5959566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize));
5969566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize));
5979566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
598676f2a66SJacob Faibussowitsch       } else {
5999566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize));
6009566063dSJacob Faibussowitsch         PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize));
6019566063dSJacob Faibussowitsch         PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
602676f2a66SJacob Faibussowitsch       }
603676f2a66SJacob Faibussowitsch     }
604676f2a66SJacob Faibussowitsch     /* Update A with merge */
605676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
606676f2a66SJacob Faibussowitsch   }
607676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
608676f2a66SJacob Faibussowitsch }
609676f2a66SJacob Faibussowitsch 
6109fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt *stacksize)
611676f2a66SJacob Faibussowitsch {
612676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
613676f2a66SJacob Faibussowitsch 
614676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
615676f2a66SJacob Faibussowitsch   while (i) {
616676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
617676f2a66SJacob Faibussowitsch 
618676f2a66SJacob Faibussowitsch     if (i == 1) {
619676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
620676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
621676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
6224d3610e3SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size, ctx) > 0) {
623676f2a66SJacob Faibussowitsch           m = stack[i].start;
624676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
6259566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*size, &l));
626676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
6279566063dSJacob 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));
628676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
6299566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size));
6309566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
631676f2a66SJacob Faibussowitsch           } else {
6329566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size));
6339566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
634676f2a66SJacob Faibussowitsch           }
635676f2a66SJacob Faibussowitsch         }
636676f2a66SJacob Faibussowitsch         /* Update A with merge */
637676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
638676f2a66SJacob Faibussowitsch         --i;
639676f2a66SJacob Faibussowitsch       }
640676f2a66SJacob Faibussowitsch     } else {
641676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
642676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
643676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
644676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
645676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
6464d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*size, arr+(stack[i-1].start)*size, ctx) > 0) {
647676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
648676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
6499566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*size, &l));
650676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
6519566063dSJacob 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));
652676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
6539566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size));
6549566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
655676f2a66SJacob Faibussowitsch             } else {
6569566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size));
6579566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
658676f2a66SJacob Faibussowitsch             }
659676f2a66SJacob Faibussowitsch           }
660676f2a66SJacob Faibussowitsch           /* Update A with merge */
661676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
662676f2a66SJacob Faibussowitsch           /* Push C up the stack */
663676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
664676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
665676f2a66SJacob Faibussowitsch         } else {
666676f2a66SJacob Faibussowitsch           /* merge C into B */
667676f2a66SJacob Faibussowitsch           mergeBC:
668676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
6694d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size, ctx) > 0) {
670676f2a66SJacob Faibussowitsch             m = stack[i].start;
671676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
6729566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-1].start, stack[i].start-1, arr+stack[i].start*size, &l));
673676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
6749566063dSJacob 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));
675676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
6769566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size));
6779566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
678676f2a66SJacob Faibussowitsch             } else {
6799566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size));
6809566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r));
681676f2a66SJacob Faibussowitsch             }
682676f2a66SJacob Faibussowitsch           }
683676f2a66SJacob Faibussowitsch           /* Update B with merge */
684676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
685676f2a66SJacob Faibussowitsch         }
686676f2a66SJacob Faibussowitsch         --i;
687676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
688676f2a66SJacob Faibussowitsch         /* merge C into B */
689676f2a66SJacob Faibussowitsch         goto mergeBC;
690676f2a66SJacob Faibussowitsch       }
691676f2a66SJacob Faibussowitsch     }
692676f2a66SJacob Faibussowitsch     if (itemp == i) break;
693676f2a66SJacob Faibussowitsch   }
694676f2a66SJacob Faibussowitsch   *stacksize = i;
695676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
696676f2a66SJacob Faibussowitsch }
697676f2a66SJacob Faibussowitsch 
6989fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortMergeCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, void *ctx, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt *stacksize)
699676f2a66SJacob Faibussowitsch {
700676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
701676f2a66SJacob Faibussowitsch 
702676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
703676f2a66SJacob Faibussowitsch   while (i) {
704676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
705676f2a66SJacob Faibussowitsch 
706676f2a66SJacob Faibussowitsch     if (i == 1) {
707676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
708676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
709676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
7104d3610e3SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize, ctx) > 0) {
711676f2a66SJacob Faibussowitsch           m = stack[i].start;
712676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
7139566063dSJacob Faibussowitsch           PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*asize, &l));
714676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
7159566063dSJacob 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));
716676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
7179566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize));
7189566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize));
7199566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
720676f2a66SJacob Faibussowitsch           } else {
7219566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize));
7229566063dSJacob Faibussowitsch             PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize));
7239566063dSJacob Faibussowitsch             PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
724676f2a66SJacob Faibussowitsch           }
725676f2a66SJacob Faibussowitsch         }
726676f2a66SJacob Faibussowitsch         /* Update A with merge */
727676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
728676f2a66SJacob Faibussowitsch         --i;
729676f2a66SJacob Faibussowitsch       }
730676f2a66SJacob Faibussowitsch     } else {
731676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
732676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
733676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
734676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
735676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
7364d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*asize, arr+(stack[i-1].start)*asize, ctx) > 0) {
737676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
738676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
7399566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*asize, &l));
740676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
7419566063dSJacob 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));
742676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
7439566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize));
7449566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize));
7459566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
746676f2a66SJacob Faibussowitsch             } else {
7479566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize));
7489566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize));
7499566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
750676f2a66SJacob Faibussowitsch             }
751676f2a66SJacob Faibussowitsch           }
752676f2a66SJacob Faibussowitsch           /* Update A with merge */
753676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
754676f2a66SJacob Faibussowitsch           /* Push C up the stack */
755676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
756676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
757676f2a66SJacob Faibussowitsch         } else {
758676f2a66SJacob Faibussowitsch           /* merge C into B */
759676f2a66SJacob Faibussowitsch           mergeBC:
760676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
7614d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize, ctx) > 0) {
762676f2a66SJacob Faibussowitsch             m = stack[i].start;
763676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
7649566063dSJacob Faibussowitsch             PetscCall(PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-1].start, stack[i].start-1, arr+stack[i].start*asize, &l));
765676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
7669566063dSJacob 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));
767676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
7689566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize));
7699566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize));
7709566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
771676f2a66SJacob Faibussowitsch             } else {
7729566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize));
7739566063dSJacob Faibussowitsch               PetscCall(PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize));
7749566063dSJacob Faibussowitsch               PetscCall(PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r));
775676f2a66SJacob Faibussowitsch             }
776676f2a66SJacob Faibussowitsch           }
777676f2a66SJacob Faibussowitsch           /* Update B with merge */
778676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
779676f2a66SJacob Faibussowitsch         }
780676f2a66SJacob Faibussowitsch         --i;
781676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
782676f2a66SJacob Faibussowitsch         /* merge C into B */
783676f2a66SJacob Faibussowitsch         goto mergeBC;
784676f2a66SJacob Faibussowitsch       }
785676f2a66SJacob Faibussowitsch     }
786676f2a66SJacob Faibussowitsch     if (itemp == i) break;
787676f2a66SJacob Faibussowitsch   }
788676f2a66SJacob Faibussowitsch   *stacksize = i;
789676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
790676f2a66SJacob Faibussowitsch }
791676f2a66SJacob Faibussowitsch 
792676f2a66SJacob Faibussowitsch /* March sequentially through the array building up a "run" of weakly increasing or strictly decreasing contiguous
793676f2a66SJacob Faibussowitsch  elements. Decreasing runs are reversed by swapping. If the run is less than minrun, artificially extend it via either
794676f2a66SJacob Faibussowitsch  binary insertion sort or regulat insertion sort */
7959fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortBuildRun_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend)
796676f2a66SJacob Faibussowitsch {
797676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
798676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
799676f2a66SJacob Faibussowitsch 
800676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
801676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
802676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
8034d3610e3SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) < 0) {
804676f2a66SJacob Faibussowitsch     ++ri;
805676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8064d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) >= 0) break;
807676f2a66SJacob Faibussowitsch       ++ri;
808676f2a66SJacob Faibussowitsch     }
809676f2a66SJacob Faibussowitsch     {
810676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
811676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
812676f2a66SJacob Faibussowitsch         COPYSWAPPY(arr+lo*size, arr+hi*size, tarr, size);
813676f2a66SJacob Faibussowitsch       }
814676f2a66SJacob Faibussowitsch     }
815676f2a66SJacob Faibussowitsch   } else {
816676f2a66SJacob Faibussowitsch     ++ri;
817676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8184d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) < 0) break;
819676f2a66SJacob Faibussowitsch       ++ri;
820676f2a66SJacob Faibussowitsch     }
821676f2a66SJacob Faibussowitsch   }
822676f2a66SJacob Faibussowitsch   if (ri < re) {
823676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
824676f2a66SJacob Faibussowitsch      binary search */
825676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
826676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8274d3610e3SJacob Faibussowitsch       PetscInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
828676f2a66SJacob Faibussowitsch     } else {
8294d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
830676f2a66SJacob Faibussowitsch     }
831676f2a66SJacob Faibussowitsch     *runend = re;
832676f2a66SJacob Faibussowitsch   } else *runend = ri;
833676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
834676f2a66SJacob Faibussowitsch }
835676f2a66SJacob Faibussowitsch 
8369fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscTimSortBuildRunWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend)
837676f2a66SJacob Faibussowitsch {
838676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
839676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
840676f2a66SJacob Faibussowitsch 
841676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
842676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
843676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
8444d3610e3SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*asize, arr+(ri*asize), ctx) < 0) {
845676f2a66SJacob Faibussowitsch     ++ri;
846676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8474d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize, ctx) >= 0) break;
848676f2a66SJacob Faibussowitsch       ++ri;
849676f2a66SJacob Faibussowitsch     }
850676f2a66SJacob Faibussowitsch     {
851676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
852676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
853676f2a66SJacob Faibussowitsch         COPYSWAPPY2(arr+lo*asize, arr+hi*asize, asize, barr+lo*bsize, barr+hi*bsize, bsize, atarr);
854676f2a66SJacob Faibussowitsch       }
855676f2a66SJacob Faibussowitsch     }
856676f2a66SJacob Faibussowitsch   } else {
857676f2a66SJacob Faibussowitsch     ++ri;
858676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8594d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize, ctx) < 0) break;
860676f2a66SJacob Faibussowitsch       ++ri;
861676f2a66SJacob Faibussowitsch     }
862676f2a66SJacob Faibussowitsch   }
863676f2a66SJacob Faibussowitsch   if (ri < re) {
864676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
865676f2a66SJacob Faibussowitsch      binary search */
866676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
867676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8684d3610e3SJacob Faibussowitsch       PetscInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
869676f2a66SJacob Faibussowitsch     } else {
8704d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
871676f2a66SJacob Faibussowitsch     }
872676f2a66SJacob Faibussowitsch     *runend = re;
873676f2a66SJacob Faibussowitsch   } else *runend = ri;
874676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
875676f2a66SJacob Faibussowitsch }
876676f2a66SJacob Faibussowitsch 
8774d3610e3SJacob Faibussowitsch /*@C
878676f2a66SJacob Faibussowitsch   PetscTimSort - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm.
879676f2a66SJacob Faibussowitsch 
880676f2a66SJacob Faibussowitsch   Not Collective
881676f2a66SJacob Faibussowitsch 
882676f2a66SJacob Faibussowitsch   Input Parameters:
883676f2a66SJacob Faibussowitsch + n    - number of values
884676f2a66SJacob Faibussowitsch . arr  - array to be sorted
885676f2a66SJacob Faibussowitsch . size - size in bytes of the datatype held in arr
8864d3610e3SJacob Faibussowitsch . cmp  - function pointer to comparison function
8874d3610e3SJacob Faibussowitsch - ctx  - optional context to be passed to comparison function, NULL if not needed
888676f2a66SJacob Faibussowitsch 
889676f2a66SJacob Faibussowitsch   Output Parameters:
890676f2a66SJacob Faibussowitsch . arr  - sorted array
891676f2a66SJacob Faibussowitsch 
8928da36234SJacob Faibussowitsch   Notes:
8938da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
8948da36234SJacob Faibussowitsch  sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims to
8958da36234SJacob Faibussowitsch  select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced arrays. To
8968da36234SJacob Faibussowitsch  do so it repeatedly triggers attempts throughout to merge adjacent runs.
8978da36234SJacob Faibussowitsch 
8988da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
8998da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
9008da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
9018da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
9028da36234SJacob Faibussowitsch   likely all contain similar data.
9038da36234SJacob Faibussowitsch 
904676f2a66SJacob Faibussowitsch   Sample usage:
905676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
906676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
907676f2a66SJacob Faibussowitsch   may also
908676f2a66SJacob Faibussowitsch  change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if
909676f2a66SJacob Faibussowitsch  the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing
9108da36234SJacob Faibussowitsch   order
911676f2a66SJacob Faibussowitsch 
912676f2a66SJacob Faibussowitsch .vb
9134d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
914676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
915676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
916676f2a66SJacob Faibussowitsch   }
917676f2a66SJacob Faibussowitsch .ve
9184d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
9194d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
920676f2a66SJacob Faibussowitsch   Then pass the function
921676f2a66SJacob Faibussowitsch .vb
9224d3610e3SJacob Faibussowitsch   PetscTimSort(n, arr, sizeof(arr[0]), my_increasing_comparison_function, ctx)
923676f2a66SJacob Faibussowitsch .ve
924676f2a66SJacob Faibussowitsch 
9254d3610e3SJacob Faibussowitsch   Fortran Notes:
9264d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
9278da36234SJacob Faibussowitsch   returns result. For example
9284d3610e3SJacob Faibussowitsch .vb
9294d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
9304d3610e3SJacob Faibussowitsch    implicit none
9314d3610e3SJacob Faibussowitsch 
9324d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
9334d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
9344d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
9354d3610e3SJacob Faibussowitsch 
9364d3610e3SJacob Faibussowitsch    if (left < right) then
9374d3610e3SJacob Faibussowitsch      result = -1
9384d3610e3SJacob Faibussowitsch    else if (left == right) then
9394d3610e3SJacob Faibussowitsch      result = 0
9404d3610e3SJacob Faibussowitsch    else
9414d3610e3SJacob Faibussowitsch      result = 1
9424d3610e3SJacob Faibussowitsch    end if
9434d3610e3SJacob Faibussowitsch    return
9444d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
9454d3610e3SJacob Faibussowitsch .ve
9464d3610e3SJacob Faibussowitsch 
947676f2a66SJacob Faibussowitsch   References:
948606c0280SSatish Balay . * - Tim Peters. https://bugs.python.org/file4451/timsort.txt
949676f2a66SJacob Faibussowitsch 
950676f2a66SJacob Faibussowitsch   Level: developer
951676f2a66SJacob Faibussowitsch 
952676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscIntSortSemiOrdered(), PetscRealSortSemiOrdered(), PetscMPIIntSortSemiOrdered()
9534d3610e3SJacob Faibussowitsch @*/
9544d3610e3SJacob Faibussowitsch PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *, void *), void *ctx)
955676f2a66SJacob Faibussowitsch {
956676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
957676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
958676f2a66SJacob Faibussowitsch   PetscTimSortBuffer buff;
959676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
960676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
961676f2a66SJacob Faibussowitsch 
962676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
963676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
964676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
965676f2a66SJacob Faibussowitsch   {
966676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
967676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
968676f2a66SJacob Faibussowitsch     while (t >= 64) {
969676f2a66SJacob Faibussowitsch       r |= t & 1;
970676f2a66SJacob Faibussowitsch       t >>= 1;
971676f2a66SJacob Faibussowitsch     }
972676f2a66SJacob Faibussowitsch     minrun = t + r;
973676f2a66SJacob Faibussowitsch   }
9744d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
9759566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun));
9764d3610e3SJacob Faibussowitsch     if (n < 64) {
9779566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "n %" PetscInt_FMT " < 64, consider using PetscSortInt() instead\n", n));
978*08401ef6SPierre Jolivet     } else PetscCheck(!(minrun < 32) && !(minrun > 65),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %" PetscInt_FMT " not in range (32,65)",minrun);
9794d3610e3SJacob Faibussowitsch   }
9809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t) minrun*size, &buff.ptr));
981676f2a66SJacob Faibussowitsch   buff.size = (size_t) minrun*size;
982676f2a66SJacob Faibussowitsch   buff.maxsize = (size_t) n*size;
983676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
984676f2a66SJacob Faibussowitsch   while (runstart < n) {
985676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
9869566063dSJacob Faibussowitsch     PetscCall(PetscTimSortBuildRun_Private((char *)arr, buff.ptr, size, cmp, ctx, n, minrun, runstart, &runend));
987676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
988676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
9899566063dSJacob Faibussowitsch     PetscCall(PetscTimSortMergeCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, &stacksize));
990676f2a66SJacob Faibussowitsch     ++stacksize;
991676f2a66SJacob Faibussowitsch     runstart = runend+1;
992676f2a66SJacob Faibussowitsch   }
993676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
994676f2a66SJacob Faibussowitsch   --stacksize;
9959566063dSJacob Faibussowitsch   PetscCall(PetscTimSortForceCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, stacksize));
9969566063dSJacob Faibussowitsch   PetscCall(PetscFree(buff.ptr));
997676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
998676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
999676f2a66SJacob Faibussowitsch }
1000676f2a66SJacob Faibussowitsch 
10014d3610e3SJacob Faibussowitsch /*@C
1002676f2a66SJacob Faibussowitsch   PetscTimSortWithArray - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm and
1003676f2a66SJacob Faibussowitsch   reorders a second array to match the first. The arrays need not be the same type.
1004676f2a66SJacob Faibussowitsch 
1005676f2a66SJacob Faibussowitsch   Not Collective
1006676f2a66SJacob Faibussowitsch 
1007676f2a66SJacob Faibussowitsch   Input Parameters:
1008676f2a66SJacob Faibussowitsch + n     - number of values
1009676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in arr
10106b867d5aSJose E. Roman . bsize - size in bytes of the datatype held in barr
10114d3610e3SJacob Faibussowitsch . cmp   - function pointer to comparison function
10124d3610e3SJacob Faibussowitsch - ctx   - optional context to be passed to comparison function, NULL if not needed
1013676f2a66SJacob Faibussowitsch 
10146b867d5aSJose E. Roman   Input/Output Parameters:
10156b867d5aSJose E. Roman + arr  - array to be sorted, on output it is sorted
10166b867d5aSJose E. Roman - barr - array to be reordered, on output it is reordered
10178da36234SJacob Faibussowitsch 
10188da36234SJacob Faibussowitsch   Notes:
10198da36234SJacob 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
10208da36234SJacob Faibussowitsch   overlap.
10218da36234SJacob Faibussowitsch 
10228da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
10238da36234SJacob Faibussowitsch   sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims
10248da36234SJacob Faibussowitsch  to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced
10258da36234SJacob Faibussowitsch  arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs.
10268da36234SJacob Faibussowitsch 
10278da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
10288da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
10298da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
10308da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
10318da36234SJacob Faibussowitsch   likely all contain similar data.
1032676f2a66SJacob Faibussowitsch 
1033676f2a66SJacob Faibussowitsch   Sample usage:
1034676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
1035676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
1036676f2a66SJacob Faibussowitsch   may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only
1037676f2a66SJacob Faibussowitsch   guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in
10388da36234SJacob Faibussowitsch   increasing order
1039676f2a66SJacob Faibussowitsch 
1040676f2a66SJacob Faibussowitsch .vb
10414d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
1042676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
1043676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
1044676f2a66SJacob Faibussowitsch   }
1045676f2a66SJacob Faibussowitsch .ve
10464d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
10474d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
1048676f2a66SJacob Faibussowitsch   Then pass the function
1049676f2a66SJacob Faibussowitsch .vb
10504d3610e3SJacob Faibussowitsch   PetscTimSortWithArray(n, arr, sizeof(arr[0]), barr, sizeof(barr[0]), my_increasing_comparison_function, ctx)
1051676f2a66SJacob Faibussowitsch .ve
1052676f2a66SJacob Faibussowitsch 
10534d3610e3SJacob Faibussowitsch   Fortran Notes:
10544d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
10558da36234SJacob Faibussowitsch   returns result. For example
10564d3610e3SJacob Faibussowitsch .vb
10574d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
10584d3610e3SJacob Faibussowitsch    implicit none
10594d3610e3SJacob Faibussowitsch 
10604d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
10614d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
10624d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
10634d3610e3SJacob Faibussowitsch 
10644d3610e3SJacob Faibussowitsch    if (left < right) then
10654d3610e3SJacob Faibussowitsch      result = -1
10664d3610e3SJacob Faibussowitsch    else if (left == right) then
10674d3610e3SJacob Faibussowitsch      result = 0
10684d3610e3SJacob Faibussowitsch    else
10694d3610e3SJacob Faibussowitsch      result = 1
10704d3610e3SJacob Faibussowitsch    end if
10714d3610e3SJacob Faibussowitsch    return
10724d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
10734d3610e3SJacob Faibussowitsch .ve
10744d3610e3SJacob Faibussowitsch 
1075676f2a66SJacob Faibussowitsch   References:
1076606c0280SSatish Balay . * - Tim Peters. https://bugs.python.org/file4451/timsort.txt
1077676f2a66SJacob Faibussowitsch 
1078676f2a66SJacob Faibussowitsch   Level: developer
1079676f2a66SJacob Faibussowitsch 
1080676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscIntSortSemiOrderedWithArray(), PetscRealSortSemiOrderedWithArrayInt(), PetscMPIIntSortSemiOrderedWithArray()
10814d3610e3SJacob Faibussowitsch @*/
10824d3610e3SJacob Faibussowitsch PetscErrorCode PetscTimSortWithArray(PetscInt n, void *arr, size_t asize, void *barr, size_t bsize, int (*cmp)(const void *, const void *, void *), void *ctx)
1083676f2a66SJacob Faibussowitsch {
1084676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
1085676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
1086676f2a66SJacob Faibussowitsch   PetscTimSortBuffer abuff, bbuff;
1087676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
1088676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
1089676f2a66SJacob Faibussowitsch 
1090676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1091676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
1092676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
1093676f2a66SJacob Faibussowitsch   {
1094676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
1095676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
1096676f2a66SJacob Faibussowitsch     while (t >= 64) {
1097676f2a66SJacob Faibussowitsch       r |= t & 1;
1098676f2a66SJacob Faibussowitsch       t >>= 1;
1099676f2a66SJacob Faibussowitsch     }
1100676f2a66SJacob Faibussowitsch     minrun = t + r;
1101676f2a66SJacob Faibussowitsch   }
11024d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
11039566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "minrun = %" PetscInt_FMT "\n", minrun));
11042c71b3e2SJacob Faibussowitsch     PetscCheckFalse((n >= 64) && ((minrun < 32) || (minrun > 65)),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %" PetscInt_FMT " not in range (32,65)",minrun);
11054d3610e3SJacob Faibussowitsch   }
11069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t) minrun*asize, &abuff.ptr));
1107676f2a66SJacob Faibussowitsch   abuff.size = (size_t) minrun*asize;
1108676f2a66SJacob Faibussowitsch   abuff.maxsize = (size_t) n*asize;
11099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((size_t) minrun*bsize, &bbuff.ptr));
1110676f2a66SJacob Faibussowitsch   bbuff.size = (size_t) minrun*bsize;
1111676f2a66SJacob Faibussowitsch   bbuff.maxsize = (size_t) n*bsize;
1112676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1113676f2a66SJacob Faibussowitsch   while (runstart < n) {
1114676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
11159566063dSJacob Faibussowitsch     PetscCall(PetscTimSortBuildRunWithArray_Private((char *)arr, abuff.ptr, asize, (char *)barr, bbuff.ptr, bsize, cmp, ctx, n, minrun, runstart, &runend));
1116676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
1117676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
11189566063dSJacob Faibussowitsch     PetscCall(PetscTimSortMergeCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, &stacksize));
1119676f2a66SJacob Faibussowitsch     ++stacksize;
1120676f2a66SJacob Faibussowitsch     runstart = runend+1;
1121676f2a66SJacob Faibussowitsch   }
1122676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
1123676f2a66SJacob Faibussowitsch   --stacksize;
11249566063dSJacob Faibussowitsch   PetscCall(PetscTimSortForceCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, stacksize));
11259566063dSJacob Faibussowitsch   PetscCall(PetscFree(abuff.ptr));
11269566063dSJacob Faibussowitsch   PetscCall(PetscFree(bbuff.ptr));
1127676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1128676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1129676f2a66SJacob Faibussowitsch }
1130676f2a66SJacob Faibussowitsch 
1131676f2a66SJacob Faibussowitsch /*@
1132676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrdered - Sorts an array of integers in place in increasing order.
1133676f2a66SJacob Faibussowitsch 
1134676f2a66SJacob Faibussowitsch    Not Collective
1135676f2a66SJacob Faibussowitsch 
1136676f2a66SJacob Faibussowitsch    Input Parameters:
1137676f2a66SJacob Faibussowitsch +  n   - number of values
1138676f2a66SJacob Faibussowitsch -  arr - array of integers
1139676f2a66SJacob Faibussowitsch 
1140676f2a66SJacob Faibussowitsch    Output Parameters:
1141676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1142676f2a66SJacob Faibussowitsch 
1143676f2a66SJacob Faibussowitsch    Notes:
1144676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortInt() is automatically used.
1145676f2a66SJacob Faibussowitsch 
1146676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortInt(). While this function works for any array of integers it is
1147676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1148a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1149676f2a66SJacob Faibussowitsch 
1150676f2a66SJacob Faibussowitsch    Level: intermediate
1151676f2a66SJacob Faibussowitsch 
1152676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortInt(), PetscSortIntWithPermutation()
1153676f2a66SJacob Faibussowitsch @*/
1154676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrdered(PetscInt n, PetscInt arr[])
1155676f2a66SJacob Faibussowitsch {
1156676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1157676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1158676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1159676f2a66SJacob Faibussowitsch   if (n < 64) {
11609566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(n, arr));
1161676f2a66SJacob Faibussowitsch   } else {
11629566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscInt), Compare_PetscInt_Private, NULL));
1163676f2a66SJacob Faibussowitsch   }
1164676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1165676f2a66SJacob Faibussowitsch }
1166676f2a66SJacob Faibussowitsch 
1167676f2a66SJacob Faibussowitsch /*@
1168676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1169676f2a66SJacob Faibussowitsch    array to match the first.
1170676f2a66SJacob Faibussowitsch 
1171676f2a66SJacob Faibussowitsch    Not Collective
1172676f2a66SJacob Faibussowitsch 
11736b867d5aSJose E. Roman    Input Parameter:
11746b867d5aSJose E. Roman .  n   - number of values
1175676f2a66SJacob Faibussowitsch 
11766b867d5aSJose E. Roman    Input/Output Parameters:
11776b867d5aSJose E. Roman +  arr1 - array of integers to be sorted, modified on output
11786b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1179676f2a66SJacob Faibussowitsch 
1180676f2a66SJacob Faibussowitsch    Notes:
1181676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1182676f2a66SJacob Faibussowitsch 
1183676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortIntWithArray(). While this function works for any array of integers it is
1184676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1185a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1186676f2a66SJacob Faibussowitsch 
1187676f2a66SJacob Faibussowitsch    Level: intermediate
1188676f2a66SJacob Faibussowitsch 
1189676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortIntWithArray(), PetscSortIntWithPermutation()
1190676f2a66SJacob Faibussowitsch @*/
1191676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt n, PetscInt arr1[], PetscInt arr2[])
1192676f2a66SJacob Faibussowitsch {
1193676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
119452e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1195676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1196676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
11976e260675SJacob Faibussowitsch   /* cannot export out to PetscIntSortWithArray here since it isn't stable */
11989566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscInt), arr2, sizeof(PetscInt), Compare_PetscInt_Private, NULL));
1199676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1200676f2a66SJacob Faibussowitsch }
1201676f2a66SJacob Faibussowitsch 
1202676f2a66SJacob Faibussowitsch /*@
1203676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrdered - Sorts an array of PetscMPIInts in place in increasing order.
1204676f2a66SJacob Faibussowitsch 
1205676f2a66SJacob Faibussowitsch    Not Collective
1206676f2a66SJacob Faibussowitsch 
1207676f2a66SJacob Faibussowitsch    Input Parameters:
1208676f2a66SJacob Faibussowitsch +  n   - number of values
1209676f2a66SJacob Faibussowitsch -  arr - array of PetscMPIInts
1210676f2a66SJacob Faibussowitsch 
1211676f2a66SJacob Faibussowitsch    Output Parameters:
1212676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1213676f2a66SJacob Faibussowitsch 
1214676f2a66SJacob Faibussowitsch    Notes:
1215676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortMPIInt() is automatically used.
1216676f2a66SJacob Faibussowitsch 
1217676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIInt(). While this function works for any array of PetscMPIInts it is
1218676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1219a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1220676f2a66SJacob Faibussowitsch 
1221676f2a66SJacob Faibussowitsch    Level: intermediate
1222676f2a66SJacob Faibussowitsch 
1223676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortMPIInt()
1224676f2a66SJacob Faibussowitsch @*/
1225676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt n, PetscMPIInt arr[])
1226676f2a66SJacob Faibussowitsch {
1227676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
122852e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1229676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1230676f2a66SJacob Faibussowitsch   if (n < 64) {
12319566063dSJacob Faibussowitsch     PetscCall(PetscSortMPIInt(n, arr));
1232676f2a66SJacob Faibussowitsch   } else {
12339566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL));
1234676f2a66SJacob Faibussowitsch   }
1235676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1236676f2a66SJacob Faibussowitsch }
1237676f2a66SJacob Faibussowitsch 
1238676f2a66SJacob Faibussowitsch /*@
1239676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1240676f2a66SJacob Faibussowitsch    array to match the first.
1241676f2a66SJacob Faibussowitsch 
1242676f2a66SJacob Faibussowitsch    Not Collective
1243676f2a66SJacob Faibussowitsch 
12446b867d5aSJose E. Roman    Input Parameter:
12456b867d5aSJose E. Roman .  n   - number of values
1246676f2a66SJacob Faibussowitsch 
12476b867d5aSJose E. Roman    Input/Output Parameters:
12486b867d5aSJose E. Roman .  arr1 - array of integers to be sorted, modified on output
12496b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1250676f2a66SJacob Faibussowitsch 
1251676f2a66SJacob Faibussowitsch    Notes:
1252676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1253676f2a66SJacob Faibussowitsch 
1254676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIIntWithArray(). While this function works for any array of integers it is
1255676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1256a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1257676f2a66SJacob Faibussowitsch 
1258676f2a66SJacob Faibussowitsch    Level: intermediate
1259676f2a66SJacob Faibussowitsch 
1260676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortMPIIntWithArray(), PetscSortMPIIntWithPermutation()
1261676f2a66SJacob Faibussowitsch @*/
1262676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt n, PetscMPIInt arr1[], PetscMPIInt arr2[])
1263676f2a66SJacob Faibussowitsch {
1264676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1265676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1266676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1267676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
12686e260675SJacob Faibussowitsch   /* cannot export out to PetscMPIIntSortWithArray here since it isn't stable */
12699566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscMPIInt), arr2, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL));
1270676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1271676f2a66SJacob Faibussowitsch }
1272676f2a66SJacob Faibussowitsch 
1273676f2a66SJacob Faibussowitsch /*@
1274676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrdered - Sorts an array of PetscReals in place in increasing order.
1275676f2a66SJacob Faibussowitsch 
1276676f2a66SJacob Faibussowitsch    Not Collective
1277676f2a66SJacob Faibussowitsch 
1278676f2a66SJacob Faibussowitsch    Input Parameters:
1279676f2a66SJacob Faibussowitsch +  n   - number of values
1280676f2a66SJacob Faibussowitsch -  arr - array of PetscReals
1281676f2a66SJacob Faibussowitsch 
1282676f2a66SJacob Faibussowitsch    Output Parameters:
1283676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1284676f2a66SJacob Faibussowitsch 
1285676f2a66SJacob Faibussowitsch    Notes:
1286676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortReal() is automatically used.
1287676f2a66SJacob Faibussowitsch 
1288676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortReal(). While this function works for any array of PetscReals it is
1289676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1290a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1291676f2a66SJacob Faibussowitsch 
1292676f2a66SJacob Faibussowitsch    Level: intermediate
1293676f2a66SJacob Faibussowitsch 
1294676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortReal(), PetscSortRealWithPermutation()
1295676f2a66SJacob Faibussowitsch @*/
1296676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrdered(PetscInt n, PetscReal arr[])
1297676f2a66SJacob Faibussowitsch {
1298676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1299676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1300676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr,2);
1301676f2a66SJacob Faibussowitsch   if (n < 64) {
13029566063dSJacob Faibussowitsch     PetscCall(PetscSortReal(n, arr));
1303676f2a66SJacob Faibussowitsch   } else {
13049566063dSJacob Faibussowitsch     PetscCall(PetscTimSort(n, arr, sizeof(PetscReal), Compare_PetscReal_Private, NULL));
1305676f2a66SJacob Faibussowitsch   }
1306676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1307676f2a66SJacob Faibussowitsch }
1308676f2a66SJacob Faibussowitsch 
1309676f2a66SJacob Faibussowitsch /*@
1310676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrderedWithArrayInt - Sorts an array of PetscReals in place in increasing order and reorders a second
1311676f2a66SJacob Faibussowitsch    array of PetscInts to match the first.
1312676f2a66SJacob Faibussowitsch 
1313676f2a66SJacob Faibussowitsch    Not Collective
1314676f2a66SJacob Faibussowitsch 
13156b867d5aSJose E. Roman    Input Parameter:
13166b867d5aSJose E. Roman .  n   - number of values
1317676f2a66SJacob Faibussowitsch 
13186b867d5aSJose E. Roman    Input/Output Parameters:
13196b867d5aSJose E. Roman .  arr1 - array of PetscReals to be sorted, modified on output
13206b867d5aSJose E. Roman -  arr2 - array of PetscReals to be reordered, modified on output
1321676f2a66SJacob Faibussowitsch 
1322676f2a66SJacob Faibussowitsch    Notes:
1323676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortRealWithArray(). While this function works for any array of PetscReals it is
1324676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1325a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1326676f2a66SJacob Faibussowitsch 
1327676f2a66SJacob Faibussowitsch    Level: intermediate
1328676f2a66SJacob Faibussowitsch 
1329676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortRealWithArrayInt(), PetscSortRealWithPermutation()
1330676f2a66SJacob Faibussowitsch @*/
1331676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt n, PetscReal arr1[], PetscInt arr2[])
1332676f2a66SJacob Faibussowitsch {
1333676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1334676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1335676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr1,2);
1336676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
13376e260675SJacob Faibussowitsch   /* cannot export out to PetscRealSortWithArrayInt here since it isn't stable */
13389566063dSJacob Faibussowitsch   PetscCall(PetscTimSortWithArray(n, arr1, sizeof(PetscReal), arr2, sizeof(PetscInt), Compare_PetscReal_Private, NULL));
1339676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1340676f2a66SJacob Faibussowitsch }
1341